OpenVDB  1.1.0
Namespaces | Classes | Enumerations | Functions
openvdb::v1_1_0::tools Namespace Reference

Namespaces

namespace  composite
 
namespace  internal
 
namespace  local
 
namespace  local_util
 
namespace  valxform
 

Classes

struct  CompReplaceOp
 
class  CsgVisitorBase
 
struct  CsgUnionVisitor
 
struct  CsgIntersectVisitor
 
struct  CsgDiffVisitor
 
class  Filter
 Filtering of VDB volumes. More...
 
struct  VectorToScalarConverter
 VectorToScalarConverter<VectorGridType>::Type is the type of a grid having the same tree configuration as VectorGridType but a scalar value type, T, where T is the type of the original vector components. More...
 
struct  ScalarToVectorConverter
 ScalarToVectorConverter<ScalarGridType>::Type is the type of a grid having the same tree configuration as ScalarGridType but value type Vec3<T> where T is ScalarGridType::ValueType. More...
 
class  Cpt
 Compute the closest-point transform of a scalar grid. More...
 
class  Curl
 Compute the curl of a scalar grid. More...
 
class  Divergence
 Computes the Divergence of a scalar grid. More...
 
class  Gradient
 Computes the Gradient of a scalar grid. More...
 
class  Laplacian
 Computes the Laplacian of a scalar grid. More...
 
class  MeanCurvature
 
class  Magnitude
 
class  Normalize
 
class  GridResampler
 
class  GridTransformer
 A GridTransformer applies a geometric transformation to an input grid using one of several sampling schemes, and stores the result in an output grid. More...
 
class  ABTransform
 This class implements the Transformer functor interface (specifically, the isAffine(), transform() and invTransform() methods) for a transform that maps an A grid into a B grid's index space such that, after resampling, A's index space and transform match B's index space and transform. More...
 
struct  PointSampler
 
struct  BoxSampler
 
struct  QuadraticSampler
 
struct  StaggeredPointSampler
 
struct  StaggeredBoxSampler
 
struct  StaggeredQuadraticSampler
 
class  GridSampler
 Base class that provides the interface for continuous sampling of values in a grid. More...
 
class  DiscreteField
 Thin wrapper class for a velocity grid. More...
 
class  EnrightField
 Analytical, divergence-free and periodic vecloity field. More...
 
class  LevelSetAdvection
 Hyperbolic advection of narrow-band level sets in an external velocity field. More...
 
class  LevelSetFilter
 Filtering (i.e. diffusion) of narrow-band level sets. More...
 
class  LevelSetFracture
 Level set fracturing. More...
 
class  LevelSetSphere
 Generates a signed distance field (or narrow band level set) to a single sphere. More...
 
class  LevelSetTracker
 Performs multi-threaded interface tracking of narrow band level sets. More...
 
class  MinMaxVoxel
 Threaded operator that finds the minimum and maximum values among the active leaf-level voxels of a grid. More...
 
class  LeafTransformer
 Threaded operator that applies a user-supplied functor to each leaf node in a LeafManager. More...
 
class  MeshToVolume
 
struct  DimToWord
 Mapping from a Log2Dim to a data type of size 2^Log2Dim bits. More...
 
struct  DimToWord< 3 >
 
struct  DimToWord< 4 >
 
struct  DimToWord< 5 >
 
struct  DimToWord< 6 >
 
class  Morphology
 
class  ParticlesToLevelSet
 
class  ParticlesToLevelSetAndId
 Use this wrapper class to convert particles into a level set and a separate index grid of closest-point particle id. The latter can be used to subsequently transfer particles attributes into separate grids. More...
 
class  ClosestPointProjector
 
class  VelocitySampler
 
class  VelocityIntegrator
 Performs runge-kutta time integration of variable order in a static velocity field. More...
 
class  PointAdvect
 
class  ConstrainedPointAdvect
 
class  UniformPointScatter
 The two point scatters UniformPointScatter and NonUniformPointScatter depend on the following two classes: More...
 
class  NonUniformPointScatter
 Non-uniform scatters of point in the active voxels. The local point count is implicitly defined as a product of of a global density and the local voxel (or tile) value. More...
 
class  PolygonPool
 Collection of quads and triangles. More...
 
class  VolumeToMesh
 Mesh any scalar grid with a continuous isosurface. More...
 

Typedefs

typedef boost::scoped_array
< openvdb::Vec3s
PointList
 Point and primitive list types. More...
 
typedef boost::scoped_array
< PolygonPool
PolygonPoolList
 Point and primitive list types. More...
 

Enumerations

enum  { GENERATE_PRIM_INDEX_GRID = 0x1 }
 Conversion flags, used to control the MeshToVolume output. More...
 
enum  { POLYFLAG_EXTERIOR = 0x1, POLYFLAG_FRACTURE_SEAM = 0x2 }
 Polygon flags, used for reference based meshing. More...
 

Functions

template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void csgUnion (GridOrTreeT &a, GridOrTreeT &b, bool prune=true)
 Given two level set grids, replace the A grid with the union of A and B. More...
 
template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void csgIntersection (GridOrTreeT &a, GridOrTreeT &b, bool prune=true)
 Given two level set grids, replace the A grid with the intersection of A and B. More...
 
template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void csgDifference (GridOrTreeT &a, GridOrTreeT &b, bool prune=true)
 Given two level set grids, replace the A grid with the difference A / B. More...
 
template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void compMax (GridOrTreeT &a, GridOrTreeT &b)
 Given grids A and B, compute max(a, b) per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty. More...
 
template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void compMin (GridOrTreeT &a, GridOrTreeT &b)
 Given grids A and B, compute min(a, b) per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty. More...
 
template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void compSum (GridOrTreeT &a, GridOrTreeT &b)
 Given grids A and B, compute a + b per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty. More...
 
template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void compMul (GridOrTreeT &a, GridOrTreeT &b)
 Given grids A and B, compute a * b per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty. More...
 
template<typename GridOrTreeT >
OPENVDB_STATIC_SPECIALIZATION void compReplace (GridOrTreeT &a, const GridOrTreeT &b)
 Copy the active voxels of B into A. More...
 
template<typename GridType , typename InterruptT >
ScalarToVectorConverter
< GridType >::Type::Ptr 
cpt (const GridType &grid, bool threaded, InterruptT *interrupt)
 Compute the Closest-Point Transform (CPT) from a distance field. More...
 
template<typename GridType >
ScalarToVectorConverter
< GridType >::Type::Ptr 
cpt (const GridType &grid, bool threaded=true)
 
template<typename GridType , typename InterruptT >
GridType::Ptr curl (const GridType &grid, bool threaded, InterruptT *interrupt)
 Compute the curl of the given vector-valued grid. More...
 
template<typename GridType >
GridType::Ptr curl (const GridType &grid, bool threaded=true)
 
template<typename GridType , typename InterruptT >
VectorToScalarConverter
< GridType >::Type::Ptr 
divergence (const GridType &grid, bool threaded, InterruptT *interrupt)
 Compute the divergence of the given vector-valued grid. More...
 
template<typename GridType >
VectorToScalarConverter
< GridType >::Type::Ptr 
divergence (const GridType &grid, bool threaded=true)
 
template<typename GridType , typename InterruptT >
ScalarToVectorConverter
< GridType >::Type::Ptr 
gradient (const GridType &grid, bool threaded, InterruptT *interrupt)
 Compute the gradient of the given scalar grid. More...
 
template<typename GridType >
ScalarToVectorConverter
< GridType >::Type::Ptr 
gradient (const GridType &grid, bool threaded=true)
 
template<typename GridType , typename InterruptT >
GridType::Ptr laplacian (const GridType &grid, bool threaded, InterruptT *interrupt)
 Compute the Laplacian of the given scalar grid. More...
 
template<typename GridType >
GridType::Ptr laplacian (const GridType &grid, bool threaded=true)
 
template<typename GridType , typename InterruptT >
GridType::Ptr meanCurvature (const GridType &grid, bool threaded, InterruptT *interrupt)
 Compute the mean curvature of the given grid. More...
 
template<typename GridType >
GridType::Ptr meanCurvature (const GridType &grid, bool threaded=true)
 
template<typename GridType , typename InterruptT >
VectorToScalarConverter
< GridType >::Type::Ptr 
magnitude (const GridType &grid, bool threaded, InterruptT *interrupt)
 Compute the magnitudes of the vectors of the given vector-valued grid. More...
 
template<typename GridType >
VectorToScalarConverter
< GridType >::Type::Ptr 
magnitude (const GridType &grid, bool threaded=true)
 
template<typename GridType , typename InterruptT >
GridType::Ptr normalize (const GridType &grid, bool threaded, InterruptT *interrupt)
 Normalize the vectors of the given vector-valued grid. More...
 
template<typename GridType >
GridType::Ptr normalize (const GridType &grid, bool threaded=true)
 
template<typename Sampler , typename Interrupter , typename GridType >
void resampleToMatch (const GridType &inGrid, GridType &outGrid, Interrupter &interrupter)
 Resample an input grid into an output grid of the same type such that, after resampling, the input and output grids coincide (apart from sampling artifacts), but the output grid's transform is unchanged. More...
 
template<typename Sampler , typename GridType >
void resampleToMatch (const GridType &inGrid, GridType &outGrid)
 Resample an input grid into an output grid of the same type such that, after resampling, the input and output grids coincide (apart from sampling artifacts), but the output grid's transform is unchanged. More...
 
template<typename Sampler , typename Interrupter , typename GridType >
void doResampleToMatch (const GridType &inGrid, GridType &outGrid, Interrupter &interrupter)
 
template<class GridType >
GridType::Ptr levelSetRebuild (const GridType &grid, float isovalue=0, float halfWidth=float(LEVEL_SET_HALF_WIDTH), const math::Transform *xform=NULL)
 Return a new grid of type GridType that contains a narrow-band level set representation of an isosurface of a given grid. More...
 
template<class GridType >
GridType::Ptr levelSetRebuild (const GridType &grid, float isovalue, float exBandWidth, float inBandWidth, const math::Transform *xform=NULL)
 Return a new grid of type GridType that contains a narrow-band level set representation of an isosurface of a given grid. More...
 
template<class GridType , typename InterruptT >
GridType::Ptr levelSetRebuild (const GridType &grid, float isovalue, float exBandWidth, float inBandWidth, const math::Transform *xform=NULL, InterruptT *interrupter=NULL)
 Return a new grid of type GridType that contains a narrow-band level set representation of an isosurface of a given grid. More...
 
template<class GridType , typename InterruptT >
boost::enable_if
< boost::is_floating_point
< typename GridType::ValueType >
, typename GridType::Ptr >
::type 
doLevelSetRebuild (const GridType &grid, typename GridType::ValueType iso, typename GridType::ValueType exWidth, typename GridType::ValueType inWidth, const math::Transform *xform, InterruptT *interrupter)
 
template<class GridType , typename InterruptT >
boost::disable_if
< boost::is_floating_point
< typename GridType::ValueType >
, typename GridType::Ptr >
::type 
doLevelSetRebuild (const GridType &, typename GridType::ValueType, typename GridType::ValueType, typename GridType::ValueType, const math::Transform *, InterruptT *)
 
template<typename GridType , typename InterruptT >
GridType::Ptr createLevelSetSphere (float radius, const openvdb::Vec3f &center, float voxelSize, float halfWidth=float(LEVEL_SET_HALF_WIDTH), InterruptT *interrupt=NULL)
 Return a grid of type GridType containing a narrow-band level set representation of a sphere. More...
 
template<typename GridType >
GridType::Ptr createLevelSetSphere (float radius, const openvdb::Vec3f &center, float voxelSize, float halfWidth=float(LEVEL_SET_HALF_WIDTH))
 Return a grid of type GridType containing a narrow-band level set representation of a sphere. More...
 
template<class GridType >
void sdfToFogVolume (GridType &grid, typename GridType::ValueType cutoffDistance=lsutilGridMax< GridType >())
 Threaded method to convert a sparse level set/SDF into a sparse fog volume. More...
 
template<class GridType >
Grid< typename
GridType::TreeType::template
ValueConverter< bool >::Type >
::Ptr 
sdfInteriorMask (const GridType &grid, typename GridType::ValueType iso=lsutilGridZero< GridType >())
 Threaded method to extract an interior region mask from a level set/SDF grid. More...
 
template<typename IterT , typename XformOp >
void foreach (const IterT &iter, XformOp &op, bool threaded=true, bool shareOp=true)
 
template<typename IterT , typename XformOp >
void foreach (const IterT &iter, const XformOp &op, bool threaded=true, bool shareOp=true)
 
template<typename InIterT , typename OutGridT , typename XformOp >
void transformValues (const InIterT &inIter, OutGridT &outGrid, XformOp &op, bool threaded=true, bool shareOp=true)
 
template<typename InIterT , typename OutGridT , typename XformOp >
void transformValues (const InIterT &inIter, OutGridT &outGrid, const XformOp &op, bool threaded=true, bool shareOp=true)
 
template<typename TreeType >
OPENVDB_STATIC_SPECIALIZATION void dilateVoxels (TreeType &tree, int count=1)
 
template<typename TreeType >
OPENVDB_STATIC_SPECIALIZATION void dilateVoxels (tree::LeafManager< TreeType > &manager, int count=1)
 
template<typename TreeType >
OPENVDB_STATIC_SPECIALIZATION void erodeVoxels (TreeType &tree, int count=1)
 
template<typename TreeType >
OPENVDB_STATIC_SPECIALIZATION void erodeVoxels (tree::LeafManager< TreeType > &manager, int count=1)
 

Typedef Documentation

typedef boost::scoped_array<openvdb::Vec3s> PointList

Point and primitive list types.

typedef boost::scoped_array<PolygonPool> PolygonPoolList

Point and primitive list types.

Enumeration Type Documentation

anonymous enum

Conversion flags, used to control the MeshToVolume output.

Enumerator
GENERATE_PRIM_INDEX_GRID 
anonymous enum

Polygon flags, used for reference based meshing.

Enumerator
POLYFLAG_EXTERIOR 
POLYFLAG_FRACTURE_SEAM 

Function Documentation

OPENVDB_STATIC_SPECIALIZATION void compMax ( GridOrTreeT &  a,
GridOrTreeT &  b 
)
inline

Given grids A and B, compute max(a, b) per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty.

OPENVDB_STATIC_SPECIALIZATION void compMin ( GridOrTreeT &  a,
GridOrTreeT &  b 
)
inline

Given grids A and B, compute min(a, b) per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty.

OPENVDB_STATIC_SPECIALIZATION void compMul ( GridOrTreeT &  a,
GridOrTreeT &  b 
)
inline

Given grids A and B, compute a * b per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty.

OPENVDB_STATIC_SPECIALIZATION void compReplace ( GridOrTreeT &  a,
const GridOrTreeT &  b 
)
inline

Copy the active voxels of B into A.

OPENVDB_STATIC_SPECIALIZATION void compSum ( GridOrTreeT &  a,
GridOrTreeT &  b 
)
inline

Given grids A and B, compute a + b per voxel (using sparse traversal). Store the result in the A grid and leave the B grid empty.

ScalarToVectorConverter< GridType >::Type::Ptr cpt ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Compute the Closest-Point Transform (CPT) from a distance field.

Returns
a new vector-valued grid with the same numerical precision as the input grid (for example, if the input grid is a DoubleGrid, the output grid will be a Vec3DGrid)
Note
The current implementation assumes all the input distance values are represented by leaf voxels and not tiles. This is true for all narrow-band level sets, which this class was originally developed for. In the future we will expand this class to also handle tile values.
ScalarToVectorConverter<GridType>::Type::Ptr openvdb::v1_1_0::tools::cpt ( const GridType &  grid,
bool  threaded = true 
)
inline
GridType::Ptr createLevelSetSphere ( float  radius,
const openvdb::Vec3f center,
float  voxelSize,
float  halfWidth = float(LEVEL_SET_HALF_WIDTH),
InterruptT *  interrupt = NULL 
)

Return a grid of type GridType containing a narrow-band level set representation of a sphere.

Parameters
radiusradius of the sphere in world units
centercenter of the sphere in world units
voxelSizevoxel size in world units
halfWidthhalf the width of the narrow band, in voxel units
interrupta pointer adhering to the util::NullInterrupter interface
Note
GridType::ValueType must be a floating-point scalar.
The leapfrog algorithm employed in this method is best suited for a single large sphere. For multiple small spheres consider using the faster algorithm in ParticlesToLevelSet.h
GridType::Ptr openvdb::v1_1_0::tools::createLevelSetSphere ( float  radius,
const openvdb::Vec3f center,
float  voxelSize,
float  halfWidth = float(LEVEL_SET_HALF_WIDTH) 
)

Return a grid of type GridType containing a narrow-band level set representation of a sphere.

Parameters
radiusradius of the sphere in world units
centercenter of the sphere in world units
voxelSizevoxel size in world units
halfWidthhalf the width of the narrow band, in voxel units
Note
GridType::ValueType must be a floating-point scalar.
The leapfrog algorithm employed in this method is best suited for a single large sphere. For multiple small spheres consider using the faster algorithm in ParticlesToLevelSet.h
OPENVDB_STATIC_SPECIALIZATION void csgDifference ( GridOrTreeT &  a,
GridOrTreeT &  b,
bool  prune = true 
)
inline

Given two level set grids, replace the A grid with the difference A / B.

Exceptions
ValueErrorif the background value of either grid is not greater than zero.
Note
This operation always leaves the B grid empty.
OPENVDB_STATIC_SPECIALIZATION void csgIntersection ( GridOrTreeT &  a,
GridOrTreeT &  b,
bool  prune = true 
)
inline

Given two level set grids, replace the A grid with the intersection of A and B.

Exceptions
ValueErrorif the background value of either grid is not greater than zero.
Note
This operation always leaves the B grid empty.
OPENVDB_STATIC_SPECIALIZATION void csgUnion ( GridOrTreeT &  a,
GridOrTreeT &  b,
bool  prune = true 
)
inline

Given two level set grids, replace the A grid with the union of A and B.

Exceptions
ValueErrorif the background value of either grid is not greater than zero.
Note
This operation always leaves the B grid empty.
GridType::Ptr curl ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Compute the curl of the given vector-valued grid.

Returns
a new vector-valued grid
GridType::Ptr openvdb::v1_1_0::tools::curl ( const GridType &  grid,
bool  threaded = true 
)
inline
OPENVDB_STATIC_SPECIALIZATION void dilateVoxels ( TreeType &  tree,
int  count = 1 
)
inline

Topologically dilate all leaf-level active voxels in the given tree, i.e., expand the set of active voxels by count voxels in the +x, -x, +y, -y, +z and -z directions, but don't change the values of any voxels, only their active states.

OPENVDB_STATIC_SPECIALIZATION void dilateVoxels ( tree::LeafManager< TreeType > &  manager,
int  count = 1 
)
inline

Topologically dilate all leaf-level active voxels in the given tree, i.e., expand the set of active voxels by count voxels in the +x, -x, +y, -y, +z and -z directions, but don't change the values of any voxels, only their active states.

VectorToScalarConverter< GridType >::Type::Ptr divergence ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Compute the divergence of the given vector-valued grid.

Returns
a new scalar-valued grid with the same numerical precision as the input grid (for example, if the input grid is a Vec3DGrid, the output grid will be a DoubleGrid)
VectorToScalarConverter<GridType>::Type::Ptr openvdb::v1_1_0::tools::divergence ( const GridType &  grid,
bool  threaded = true 
)
inline
boost::enable_if<boost::is_floating_point<typename GridType::ValueType>,typename GridType::Ptr>::type openvdb::v1_1_0::tools::doLevelSetRebuild ( const GridType &  grid,
typename GridType::ValueType  iso,
typename GridType::ValueType  exWidth,
typename GridType::ValueType  inWidth,
const math::Transform *  xform,
InterruptT *  interrupter 
)
inline

The normal entry points for level set rebuild are the levelSetRebuild() functions. doLevelSetRebuild() is mainly for internal use, but when the isovalue and half band widths are given in ValueType units (for example, if they are queried from a grid), it might be more convenient to call this function directly.

boost::disable_if<boost::is_floating_point<typename GridType::ValueType>,typename GridType::Ptr>::type openvdb::v1_1_0::tools::doLevelSetRebuild ( const GridType &  ,
typename GridType::ValueType  ,
typename GridType::ValueType  ,
typename GridType::ValueType  ,
const math::Transform *  ,
InterruptT *   
)
inline
void openvdb::v1_1_0::tools::doResampleToMatch ( const GridType &  inGrid,
GridType &  outGrid,
Interrupter &  interrupter 
)
inline

The normal entry points for resampling are the resampleToMatch() functions, which correctly handle level set grids under scaling and shearing. doResampleToMatch() is mainly for internal use but is typically faster for level sets, and correct provided that no scaling or shearing is needed.

Warning
Do not use this function to scale or shear a level set grid.
OPENVDB_STATIC_SPECIALIZATION void erodeVoxels ( TreeType &  tree,
int  count = 1 
)
inline

Topologically erode all leaf-level active voxels in the given tree, i.e., shrink the set of active voxels by count voxels in the +x, -x, +y, -y, +z and -z directions, but don't change the values of any voxels, only their active states.

OPENVDB_STATIC_SPECIALIZATION void erodeVoxels ( tree::LeafManager< TreeType > &  manager,
int  count = 1 
)
inline

Topologically erode all leaf-level active voxels in the given tree, i.e., shrink the set of active voxels by count voxels in the +x, -x, +y, -y, +z and -z directions, but don't change the values of any voxels, only their active states.

void foreach ( const IterT &  iter,
XformOp &  op,
bool  threaded = true,
bool  shareOp = true 
)
inline

Iterate over a grid and at each step call op(iter).

Parameters
iteran iterator over a grid or its tree (Grid::ValueOnCIter, Tree::NodeIter, etc.)
opa functor of the form void op(const IterT&), where IterT is the type of iter
threadedif true, transform multiple values of the grid in parallel
shareOpif true and threaded is true, all threads use the same functor; otherwise, each thread gets its own copy of the functor
Example:
Multiply all values (both set and unset) of a scalar, floating-point grid by two.
struct Local {
static inline void op(const FloatGrid::ValueAllIter& iter) {
iter.setValue(*iter * 2);
}
};
FloatGrid grid;
tools::foreach(grid.beginValueAll(), Local::op);
Example:
Rotate all active vectors of a vector grid by 45 degrees about the y axis.
namespace {
struct MatMul {
MatMul(const math::Mat3s& mat): M(mat) {}
inline void operator()(const VectorGrid::ValueOnIter& iter) const {
iter.setValue(M.transform(*iter));
}
};
}
{
VectorGrid grid;
tools::foreach(grid.beginValueOn(),
MatMul(math::rotation<math::Mat3s>(math::Y, M_PI_4)));
}
Note
For more complex operations that require finer control over threading, consider using tbb::parallel_for() or tbb::parallel_reduce() in conjunction with a tree::IteratorRange that wraps a grid or tree iterator.
void foreach ( const IterT &  iter,
const XformOp &  op,
bool  threaded = true,
bool  shareOp = true 
)
inline
ScalarToVectorConverter< GridType >::Type::Ptr gradient ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Compute the gradient of the given scalar grid.

Returns
a new vector-valued grid with the same numerical precision as the input grid (for example, if the input grid is a DoubleGrid, the output grid will be a Vec3DGrid)
ScalarToVectorConverter<GridType>::Type::Ptr openvdb::v1_1_0::tools::gradient ( const GridType &  grid,
bool  threaded = true 
)
inline
GridType::Ptr laplacian ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Compute the Laplacian of the given scalar grid.

Returns
a new scalar grid
GridType::Ptr openvdb::v1_1_0::tools::laplacian ( const GridType &  grid,
bool  threaded = true 
)
inline
GridType::Ptr levelSetRebuild ( const GridType &  grid,
float  isovalue = 0,
float  halfWidth = float(LEVEL_SET_HALF_WIDTH),
const math::Transform *  xform = NULL 
)
inline

Return a new grid of type GridType that contains a narrow-band level set representation of an isosurface of a given grid.

Parameters
grida scalar, floating-point grid with one or more disjoint, closed isosurfaces at the given isovalue
isovaluethe isovalue that defines the implicit surface (defaults to zero, which is typical if the input grid is already a level set or a SDF).
halfWidthhalf the width of the narrow band, in voxel units (defaults to 3 voxels, which is required for some level set operations)
xformoptional transform for the output grid (if not provided, the transform of the input grid will be matched)
Exceptions
TypeErrorif grid is not scalar or not floating-point
Note
If the input grid contains overlapping isosurfaces, interior edges will be lost.
GridType::Ptr levelSetRebuild ( const GridType &  grid,
float  isovalue,
float  exBandWidth,
float  inBandWidth,
const math::Transform *  xform = NULL 
)
inline

Return a new grid of type GridType that contains a narrow-band level set representation of an isosurface of a given grid.

Parameters
grida scalar, floating-point grid with one or more disjoint, closed isosurfaces at the given isovalue
isovaluethe isovalue that defines the implicit surface
exBandWidththe exterior narrow-band width in voxel units
inBandWidththe interior narrow-band width in voxel units
xformoptional transform for the output grid (if not provided, the transform of the input grid will be matched)
Exceptions
TypeErrorif grid is not scalar or not floating-point
Note
If the input grid contains overlapping isosurfaces, interior edges will be lost.
GridType::Ptr levelSetRebuild ( const GridType &  grid,
float  isovalue,
float  exBandWidth,
float  inBandWidth,
const math::Transform *  xform = NULL,
InterruptT *  interrupter = NULL 
)
inline

Return a new grid of type GridType that contains a narrow-band level set representation of an isosurface of a given grid.

Parameters
grida scalar, floating-point grid with one or more disjoint, closed isosurfaces at the given isovalue
isovaluethe isovalue that defines the implicit surface
exBandWidththe exterior narrow-band width in voxel units
inBandWidththe interior narrow-band width in voxel units
xformoptional transform for the output grid (if not provided, the transform of the input grid will be matched)
interrupteroptional interrupter object
Exceptions
TypeErrorif grid is not scalar or not floating-point
Note
If the input grid contains overlapping isosurfaces, interior edges will be lost.
VectorToScalarConverter< GridType >::Type::Ptr magnitude ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Compute the magnitudes of the vectors of the given vector-valued grid.

Returns
a new scalar-valued grid with the same numerical precision as the input grid (for example, if the input grid is a Vec3DGrid, the output grid will be a DoubleGrid)
VectorToScalarConverter<GridType>::Type::Ptr openvdb::v1_1_0::tools::magnitude ( const GridType &  grid,
bool  threaded = true 
)
inline
GridType::Ptr meanCurvature ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Compute the mean curvature of the given grid.

Returns
a new grid
GridType::Ptr openvdb::v1_1_0::tools::meanCurvature ( const GridType &  grid,
bool  threaded = true 
)
inline
GridType::Ptr normalize ( const GridType &  grid,
bool  threaded,
InterruptT *  interrupt 
)
inline

Normalize the vectors of the given vector-valued grid.

Returns
a new vector-valued grid
GridType::Ptr openvdb::v1_1_0::tools::normalize ( const GridType &  grid,
bool  threaded = true 
)
inline
void resampleToMatch ( const GridType &  inGrid,
GridType &  outGrid,
Interrupter &  interrupter 
)
inline

Resample an input grid into an output grid of the same type such that, after resampling, the input and output grids coincide (apart from sampling artifacts), but the output grid's transform is unchanged.

Specifically, this function resamples the input grid into the output grid's index space, using a sampling kernel like PointSampler, BoxSampler, or QuadraticSampler.

Parameters
inGridthe grid to be resampled
outGridthe grid into which to write the resampled voxel data
interrupteran object adhering to the util::NullInterrupter interface
Example:
// Create an input grid with the default identity transform
// and populate it with a level-set sphere.
FloatGrid::ConstPtr src = tools::makeSphere(...);
// Create an output grid and give it a uniform-scale transform.
const float voxelSize = 0.5;
dest->setTransform(math::Transform::createLinearTransform(voxelSize));
// Resample the input grid into the output grid, reproducing
// the level-set sphere at a smaller voxel size.
MyInterrupter interrupter = ...;
tools::resampleToMatch<tools::QuadraticSampler>(*src, *dest, interrupter);
void resampleToMatch ( const GridType &  inGrid,
GridType &  outGrid 
)
inline

Resample an input grid into an output grid of the same type such that, after resampling, the input and output grids coincide (apart from sampling artifacts), but the output grid's transform is unchanged.

Specifically, this function resamples the input grid into the output grid's index space, using a sampling kernel like PointSampler, BoxSampler, or QuadraticSampler.

Parameters
inGridthe grid to be resampled
outGridthe grid into which to write the resampled voxel data
Example:
// Create an input grid with the default identity transform
// and populate it with a level-set sphere.
FloatGrid::ConstPtr src = tools::makeSphere(...);
// Create an output grid and give it a uniform-scale transform.
const float voxelSize = 0.5;
dest->setTransform(math::Transform::createLinearTransform(voxelSize));
// Resample the input grid into the output grid, reproducing
// the level-set sphere at a smaller voxel size.
tools::resampleToMatch<tools::QuadraticSampler>(*src, *dest);
Grid< typename GridType::TreeType::template ValueConverter< bool >::Type >::Ptr sdfInteriorMask ( const GridType &  grid,
typename GridType::ValueType  iso = lsutilGridZero<GridType>() 
)
inline

Threaded method to extract an interior region mask from a level set/SDF grid.

Returns
a shared pointer to a new boolean grid with the same tree configuration and transform as the incoming grid and whose active voxels correspond to the interior of the input SDF
Parameters
grida level set/SDF grid
isothreshold below which values are considered to be part of the interior region
void sdfToFogVolume ( GridType &  grid,
typename GridType::ValueType  cutoffDistance = lsutilGridMax<GridType>() 
)
inline

Threaded method to convert a sparse level set/SDF into a sparse fog volume.

For a level set, the active and negative-valued interior half of the narrow band becomes a linear ramp from 0 to 1; the inactive interior becomes active with a constant value of 1; and the exterior, including the background and the active exterior half of the narrow band, becomes inactive with a constant value of 0. The interior, though active, remains sparse.

For a generic SDF, a specified cutoff distance determines the width of the ramp, but otherwise the result is the same as for a level set.

Parameters
gridlevel set/SDF grid to transform
cutoffDistanceoptional world space cutoff distance for the ramp (automatically clamped if greater than the interior narrow band width)
void transformValues ( const InIterT &  inIter,
OutGridT &  outGrid,
XformOp &  op,
bool  threaded = true,
bool  shareOp = true 
)
inline

Iterate over a grid and at each step call op(iter, accessor) to populate (via the accessor) the given output grid, whose ValueType need not be the same as the input grid's.

Parameters
inItera non-const or (preferably) const iterator over an input grid or its tree (Grid::ValueOnCIter, Tree::NodeIter, etc.)
outGridan empty grid to be populated
opa functor of the form void op(const InIterT&, OutGridT::ValueAccessor&), where InIterT is the type of inIter
threadedif true, transform multiple values of the input grid in parallel
shareOpif true and threaded is true, all threads use the same functor; otherwise, each thread gets its own copy of the functor
Example:
Populate a scalar floating-point grid with the lengths of the vectors from all active voxels of a vector-valued input grid.
struct Local {
static void op(
FloatGrid::ValueAccessor& accessor)
{
if (iter.isVoxelValue()) { // set a single voxel
accessor.setValue(iter.getCoord(), iter->length());
} else { // fill an entire tile
CoordBBox bbox;
iter.getBoundingBox(bbox);
accessor.getTree()->fill(bbox, iter->length());
}
}
};
Vec3fGrid inGrid;
FloatGrid outGrid;
tools::transformValues(inGrid.cbeginValueOn(), outGrid, Local::op);
Note
For more complex operations that require finer control over threading, consider using tbb::parallel_for() or tbb::parallel_reduce() in conjunction with a tree::IteratorRange that wraps a grid or tree iterator.
void transformValues ( const InIterT &  inIter,
OutGridT &  outGrid,
const XformOp &  op,
bool  threaded = true,
bool  shareOp = true 
)
inline