OpenVDB  11.0.0
Classes | Namespaces | Typedefs
PoissonSolver.h File Reference

Solve Poisson's equation ∇2x = b for x, where b is a vector comprising the values of all of the active voxels in a grid. More...

#include <openvdb/Types.h>
#include <openvdb/math/ConjGradient.h>
#include <openvdb/tree/LeafManager.h>
#include <openvdb/tree/Tree.h>
#include <openvdb/util/NullInterrupter.h>
#include "Morphology.h"
#include <openvdb/openvdb.h>

Go to the source code of this file.

Classes

struct  DirichletBoundaryOp< ValueType >
 Dirichlet boundary condition functor. More...
 

Namespaces

 openvdb
 
 openvdb::v11_0
 
 openvdb::v11_0::tools
 
 openvdb::v11_0::tools::poisson
 

Typedefs

using VIndex = Int32
 
using LaplacianMatrix = math::pcg::SparseStencilMatrix< double, 7 >
 The type of a matrix used to represent a three-dimensional Laplacian operator. More...
 

Functions

template<typename TreeType >
TreeType::Ptr solve (const TreeType &, math::pcg::State &, bool staggered=false)
 Solve ∇2x = b for x, where b is a vector comprising the values of all of the active voxels in the input tree. More...
 
template<typename TreeType , typename Interrupter >
TreeType::Ptr solve (const TreeType &, math::pcg::State &, Interrupter &, bool staggered=false)
 Solve ∇2x = b for x, where b is a vector comprising the values of all of the active voxels in the input tree. More...
 
template<typename TreeType , typename BoundaryOp , typename Interrupter >
TreeType::Ptr solveWithBoundaryConditions (const TreeType &, const BoundaryOp &, math::pcg::State &, Interrupter &, bool staggered=false)
 Solve ∇2x = b for x with user-specified boundary conditions, where b is a vector comprising the values of all of the active voxels in the input tree or domain mask if provided. More...
 
template<typename PreconditionerType , typename TreeType , typename BoundaryOp , typename Interrupter >
TreeType::Ptr solveWithBoundaryConditionsAndPreconditioner (const TreeType &, const BoundaryOp &, math::pcg::State &, Interrupter &, bool staggered=false)
 Solve ∇2x = b for x with user-specified boundary conditions, where b is a vector comprising the values of all of the active voxels in the input tree or domain mask if provided. More...
 
template<typename PreconditionerType , typename TreeType , typename DomainTreeType , typename BoundaryOp , typename Interrupter >
TreeType::Ptr solveWithBoundaryConditionsAndPreconditioner (const TreeType &, const DomainTreeType &, const BoundaryOp &, math::pcg::State &, Interrupter &, bool staggered=false)
 Solve ∇2x = b for x with user-specified boundary conditions, where b is a vector comprising the values of all of the active voxels in the input tree or domain mask if provided. More...
 
Low-level functions
template<typename VIndexTreeType >
void populateIndexTree (VIndexTreeType &)
 Overwrite each active voxel in the given scalar tree with a sequential index, starting from zero. More...
 
template<typename TreeType >
TreeType::template ValueConverter< VIndex >::Type::Ptr createIndexTree (const TreeType &)
 Iterate over the active voxels of the input tree and for each one assign its index in the iteration sequence to the corresponding voxel of an integer-valued output tree. More...
 
template<typename VectorValueType , typename SourceTreeType >
math::pcg::Vector< VectorValueType >::Ptr createVectorFromTree (const SourceTreeType &source, const typename SourceTreeType::template ValueConverter< VIndex >::Type &index)
 Return a vector of the active voxel values of the scalar-valued source tree. More...
 
template<typename TreeValueType , typename VIndexTreeType , typename VectorValueType >
VIndexTreeType::template ValueConverter< TreeValueType >::Type::Ptr createTreeFromVector (const math::pcg::Vector< VectorValueType > &values, const VIndexTreeType &index, const TreeValueType &background)
 Return a tree with the same active voxel topology as the index tree but whose voxel values are taken from the the given vector. More...
 
template<typename BoolTreeType >
LaplacianMatrix::Ptr createISLaplacian (const typename BoolTreeType::template ValueConverter< VIndex >::Type &vectorIndexTree, const BoolTreeType &interiorMask, bool staggered=false)
 Generate a sparse matrix of the index-space (Δx = 1) Laplacian operator using second-order finite differences. More...
 
template<typename BoolTreeType , typename BoundaryOp >
LaplacianMatrix::Ptr createISLaplacianWithBoundaryConditions (const typename BoolTreeType::template ValueConverter< VIndex >::Type &vectorIndexTree, const BoolTreeType &interiorMask, const BoundaryOp &boundaryOp, typename math::pcg::Vector< LaplacianMatrix::ValueType > &source, bool staggered=false)
 Generate a sparse matrix of the index-space (Δx = 1) Laplacian operator with user-specified boundary conditions using second-order finite differences. More...
 

Detailed Description

Solve Poisson's equation ∇2x = b for x, where b is a vector comprising the values of all of the active voxels in a grid.

Authors
D.J. Hill, Peter Cucka
Example:
Solve for the pressure in a cubic tank of liquid, assuming uniform boundary conditions:
FloatTree source(/*background=*/0.0f);
// Activate voxels to indicate that they contain liquid.
source.fill(CoordBBox(Coord(0, -10, 0), Coord(10, 0, 10)), /*value=*/0.0f);
math::pcg::State state = math::pcg::terminationDefaults<float>();
FloatTree::Ptr solution = tools::poisson::solve(source, state);
Example:
Solve for the pressure, P, in a cubic tank of liquid that is open at the top. Boundary conditions are P = 0 at the top, ∂P/∂y = −1 at the bottom and ∂P/∂x = 0 at the sides:
               P = 0
            +--------+ (N,0,N)
           /|       /|
  (0,0,0) +--------+ |
          | |      | | dP/dx = 0
dP/dx = 0 | +------|-+
          |/       |/
 (0,-N,0) +--------+ (N,-N,0)
          dP/dy = -1
const int N = 10;
DoubleTree source(/*background=*/0.0);
// Activate voxels to indicate that they contain liquid.
source.fill(CoordBBox(Coord(0, -N, 0), Coord(N, 0, N)), /*value=*/0.0);
auto boundary = [](const openvdb::Coord& ijk, const openvdb::Coord& neighbor,
double& source, double& diagonal)
{
if (neighbor.x() == ijk.x() && neighbor.z() == ijk.z()) {
if (neighbor.y() < ijk.y()) source -= 1.0;
else diagonal -= 1.0;
}
};
math::pcg::State state = math::pcg::terminationDefaults<double>();
util::NullInterrupter interrupter;
source, boundary, state, interrupter);