OpenVDB  12.0.0
Classes | Typedefs
openvdb::v12_0::tools::poisson Namespace Reference

Classes

struct  DirichletBoundaryOp
 Dirichlet boundary condition functor. More...
 

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...
 

Typedef Documentation

The type of a matrix used to represent a three-dimensional Laplacian operator.

using VIndex = Int32

Function Documentation

TreeType::template ValueConverter< VIndex >::Type::Ptr createIndexTree ( const TreeType &  tree)
inline

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.

LaplacianMatrix::Ptr createISLaplacian ( const typename BoolTreeType::template ValueConverter< VIndex >::Type &  vectorIndexTree,
const BoolTreeType &  interiorMask,
bool  staggered = false 
)
inline

Generate a sparse matrix of the index-space (Δx = 1) Laplacian operator using second-order finite differences.

This construction assumes homogeneous Dirichlet boundary conditions (exterior grid points are zero).

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 
)
inline

Generate a sparse matrix of the index-space (Δx = 1) Laplacian operator with user-specified boundary conditions using second-order finite differences.

Each thread gets its own copy of boundaryOp, which should be a functor of the form

struct BoundaryOp {
using ValueType = LaplacianMatrix::ValueType;
void operator()(
const Coord& ijk, // coordinates of a boundary voxel
const Coord& ijkNeighbor, // coordinates of an exterior neighbor of ijk
ValueType& source, // element of source vector corresponding to ijk
ValueType& diagonal // element of Laplacian matrix corresponding to ijk
) const;
};

The functor is called for each of the exterior neighbors of each boundary voxel (i,&nbsp j,&nbsp k), and it must specify a boundary condition for (i,&nbsp j,&nbsp k) by modifying one or both of two provided values: an entry in the given source vector corresponding to (i,&nbsp j,&nbsp k) and the weighting coefficient for (i,&nbsp j,&nbsp k) in the Laplacian matrix.

VIndexTreeType::template ValueConverter< TreeValueType >::Type::Ptr createTreeFromVector ( const math::pcg::Vector< VectorValueType > &  values,
const VIndexTreeType &  index,
const TreeValueType &  background 
)
inline

Return a tree with the same active voxel topology as the index tree but whose voxel values are taken from the the given vector.

The voxel whose value in the index tree is n gets assigned the nth element of the vector.

Parameters
indexa tree with value type VIndex that maps voxels to elements of values
valuesa vector of values with which to populate the active voxels of the output tree
backgroundthe value for the inactive voxels of the output tree
math::pcg::Vector< VectorValueType >::Ptr createVectorFromTree ( const SourceTreeType &  source,
const typename SourceTreeType::template ValueConverter< VIndex >::Type &  index 
)
inline

Return a vector of the active voxel values of the scalar-valued source tree.

The nth element of the vector corresponds to the voxel whose value in the index tree is n.

Parameters
sourcea tree with a scalar value type
indexa tree of the same configuration as source but with value type VIndex that maps voxels to elements of the output vector
void populateIndexTree ( VIndexTreeType &  result)
inline

Overwrite each active voxel in the given scalar tree with a sequential index, starting from zero.

TreeType::Ptr solve ( const TreeType &  inTree,
math::pcg::State 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.

Returns
a new tree, with the same active voxel topology as the input tree, whose voxel values are the elements of the solution vector x.

On input, the State object should specify convergence criteria (minimum error and maximum number of iterations); on output, it gives the actual termination conditions.

The solution is computed using the conjugate gradient method with (where possible) incomplete Cholesky preconditioning, falling back to Jacobi preconditioning.

See also
solveWithBoundaryConditions
TreeType::Ptr solve ( const TreeType &  inTree,
math::pcg::State state,
Interrupter &  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.

Returns
a new tree, with the same active voxel topology as the input tree, whose voxel values are the elements of the solution vector x.

On input, the State object should specify convergence criteria (minimum error and maximum number of iterations); on output, it gives the actual termination conditions.

The solution is computed using the conjugate gradient method with (where possible) incomplete Cholesky preconditioning, falling back to Jacobi preconditioning.

See also
solveWithBoundaryConditions
TreeType::Ptr solveWithBoundaryConditions ( const TreeType &  inTree,
const BoundaryOp &  boundaryOp,
math::pcg::State state,
Interrupter &  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.

Returns
a new tree, with the same active voxel topology as the input tree, whose voxel values are the elements of the solution vector x.

On input, the State object should specify convergence criteria (minimum error and maximum number of iterations); on output, it gives the actual termination conditions.

The solution is computed using the conjugate gradient method with the specified type of preconditioner (default: incomplete Cholesky), falling back to Jacobi preconditioning if necessary.

Each thread gets its own copy of the BoundaryOp, which should be a functor of the form

struct BoundaryOp {
using ValueType = LaplacianMatrix::ValueType;
void operator()(
const Coord& ijk, // coordinates of a boundary voxel
const Coord& ijkNeighbor, // coordinates of an exterior neighbor of ijk
ValueType& source, // element of b corresponding to ijk
ValueType& diagonal // element of Laplacian matrix corresponding to ijk
) const;
};

The functor is called for each of the exterior neighbors of each boundary voxel (i,&nbsp j,&nbsp k), and it must specify a boundary condition for (i,&nbsp j,&nbsp k) by modifying one or both of two provided values: the entry in the source vector b corresponding to (i,&nbsp j,&nbsp k) and the weighting coefficient for (i,&nbsp j,&nbsp k) in the Laplacian operator matrix.

See also
solve
TreeType::Ptr solveWithBoundaryConditionsAndPreconditioner ( const TreeType &  inTree,
const BoundaryOp &  boundaryOp,
math::pcg::State state,
Interrupter &  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.

Returns
a new tree, with the same active voxel topology as the input tree, whose voxel values are the elements of the solution vector x.

On input, the State object should specify convergence criteria (minimum error and maximum number of iterations); on output, it gives the actual termination conditions.

The solution is computed using the conjugate gradient method with the specified type of preconditioner (default: incomplete Cholesky), falling back to Jacobi preconditioning if necessary.

Each thread gets its own copy of the BoundaryOp, which should be a functor of the form

struct BoundaryOp {
using ValueType = LaplacianMatrix::ValueType;
void operator()(
const Coord& ijk, // coordinates of a boundary voxel
const Coord& ijkNeighbor, // coordinates of an exterior neighbor of ijk
ValueType& source, // element of b corresponding to ijk
ValueType& diagonal // element of Laplacian matrix corresponding to ijk
) const;
};

The functor is called for each of the exterior neighbors of each boundary voxel (i,&nbsp j,&nbsp k), and it must specify a boundary condition for (i,&nbsp j,&nbsp k) by modifying one or both of two provided values: the entry in the source vector b corresponding to (i,&nbsp j,&nbsp k) and the weighting coefficient for (i,&nbsp j,&nbsp k) in the Laplacian operator matrix.

See also
solve
TreeType::Ptr solveWithBoundaryConditionsAndPreconditioner ( const TreeType &  inTree,
const DomainTreeType &  domainMask,
const BoundaryOp &  boundaryOp,
math::pcg::State state,
Interrupter &  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.

Returns
a new tree, with the same active voxel topology as the input tree, whose voxel values are the elements of the solution vector x.

On input, the State object should specify convergence criteria (minimum error and maximum number of iterations); on output, it gives the actual termination conditions.

The solution is computed using the conjugate gradient method with the specified type of preconditioner (default: incomplete Cholesky), falling back to Jacobi preconditioning if necessary.

Each thread gets its own copy of the BoundaryOp, which should be a functor of the form

struct BoundaryOp {
using ValueType = LaplacianMatrix::ValueType;
void operator()(
const Coord& ijk, // coordinates of a boundary voxel
const Coord& ijkNeighbor, // coordinates of an exterior neighbor of ijk
ValueType& source, // element of b corresponding to ijk
ValueType& diagonal // element of Laplacian matrix corresponding to ijk
) const;
};

The functor is called for each of the exterior neighbors of each boundary voxel (i,&nbsp j,&nbsp k), and it must specify a boundary condition for (i,&nbsp j,&nbsp k) by modifying one or both of two provided values: the entry in the source vector b corresponding to (i,&nbsp j,&nbsp k) and the weighting coefficient for (i,&nbsp j,&nbsp k) in the Laplacian operator matrix.

See also
solve