OpenVDB  9.0.1
Classes | Public Member Functions | List of all members
VolumeAdvection< VelocityGridT, StaggeredVelocity, InterrupterType > Class Template Reference

Performs advections of an arbitrary type of volume in a static velocity field. The advections are performed by means of various derivatives of Semi-Lagrangian integration, i.e. backwards tracking along the hyperbolic characteristics followed by interpolation. More...

#include <openvdb/tools/VolumeAdvect.h>

Public Member Functions

 VolumeAdvection (const VelocityGridT &velGrid, InterrupterType *interrupter=nullptr)
 Constructor. More...
 
virtual ~VolumeAdvection ()
 
int spatialOrder () const
 Return the spatial order of accuracy of the advection scheme. More...
 
int temporalOrder () const
 Return the temporal order of accuracy of the advection scheme. More...
 
void setIntegrator (Scheme::SemiLagrangian integrator)
 Set the integrator (see details in the table above) More...
 
Scheme::SemiLagrangian getIntegrator () const
 Return the integrator (see details in the table above) More...
 
void setLimiter (Scheme::Limiter limiter)
 Set the limiter (see details above) More...
 
Scheme::Limiter getLimiter () const
 Retrun the limiter (see details above) More...
 
bool isLimiterOn () const
 Return true if a limiter will be applied based on the current settings. More...
 
size_t getGrainSize () const
 
void setGrainSize (size_t grainsize)
 Set the grain-size used for multi-threading. More...
 
int getSubSteps () const
 
void setSubSteps (int substeps)
 Set the number of sub-steps per integration. More...
 
double getMaxVelocity () const
 Return the maximum magnitude of the velocity in the advection velocity field defined during construction. More...
 
template<typename VolumeGridT >
int getMaxDistance (const VolumeGridT &inGrid, double dt) const
 
template<typename VolumeGridT , typename VolumeSamplerT >
VolumeGridT::Ptr advect (const VolumeGridT &inGrid, double timeStep)
 
template<typename VolumeGridT , typename MaskGridT , typename VolumeSamplerT >
VolumeGridT::Ptr advect (const VolumeGridT &inGrid, const MaskGridT &mask, double timeStep)
 

Detailed Description

template<typename VelocityGridT = Vec3fGrid, bool StaggeredVelocity = false, typename InterrupterType = util::NullInterrupter>
class openvdb::v9_0::tools::VolumeAdvection< VelocityGridT, StaggeredVelocity, InterrupterType >

Performs advections of an arbitrary type of volume in a static velocity field. The advections are performed by means of various derivatives of Semi-Lagrangian integration, i.e. backwards tracking along the hyperbolic characteristics followed by interpolation.

Note
Optionally a limiter can be combined with the higher-order integration schemes MacCormack and BFECC. There are two types of limiters (CLAMP and REVERT) that suppress non-physical oscillations by means of either claminging or reverting to a first-order schemes when the function is not bounded by the cell values used for tri-linear interpolation.
The supported integrations schemes:
///
///    ================================================================
///    |  Lable | Accuracy |  Integration Scheme   |  Interpolations  |
///    |        |Time/Space|                       |  velocity/volume |
///    ================================================================
///    |  SEMI  |   1/1    | Semi-Lagrangian       |        1/1       |
///    |  MID   |   2/1    | Mid-Point             |        2/1       |
///    |  RK3   |   3/1    | 3rd Order Runge-Kutta |        3/1       |
///    |  RK4   |   4/1    | 4th Order Runge-Kutta |        4/1       |
///    |  MAC   |   2/2    | MacCormack            |        2/2       |
///    |  BFECC |   2/2    | BFECC                 |        3/2       |
///    ================================================================
/// 

Constructor & Destructor Documentation

VolumeAdvection ( const VelocityGridT &  velGrid,
InterrupterType *  interrupter = nullptr 
)
inline

Constructor.

Parameters
velGridVelocity grid responsible for the (passive) advection.
interrupterOptional interrupter used to prematurely end computations.
Note
The velocity field is assumed to be constant for the duration of the advection.
virtual ~VolumeAdvection ( )
inlinevirtual

Member Function Documentation

VolumeGridT::Ptr advect ( const VolumeGridT &  inGrid,
double  timeStep 
)
inline
Returns
Returns a new grid that is the result of passive advection of all the active values the input grid by timeStep.
Parameters
inGridThe input grid to be advected (unmodified)
timeStepTime-step of the Runge-Kutta integrator.

This method will advect all of the active values in the input inGrid. To achieve this a deep-copy is dilated to account for the material transport. This dilation step can be slow for large time steps dt or a velocity field with large magnitudes.

Warning
If the VolumeSamplerT is of higher order than one (i.e. tri-linear interpolation) instabilities are known to occur. To suppress those monotonicity constrains or flux-limiters need to be applies.
Exceptions
RuntimeErrorif inGrid does not have uniform voxels.
VolumeGridT::Ptr advect ( const VolumeGridT &  inGrid,
const MaskGridT &  mask,
double  timeStep 
)
inline
Returns
Returns a new grid that is the result of passive advection of the active values in inGrid that intersect the active values in mask. The time of the output grid is incremented by timeStep.
Parameters
inGridThe input grid to be advected (unmodified).
maskThe mask of active values defining the active voxels in inGrid on which to perform advection. Only if a value is active in both grids will it be modified.
timeStepTime-step for a single Runge-Kutta integration step.

This method will advect all of the active values in the input inGrid that intersects with the active values in mask. To achieve this a deep-copy is dilated to account for the material transport and finally cropped to the intersection with mask. The dilation step can be slow for large time steps dt or fast moving velocity fields.

Warning
If the VolumeSamplerT is of higher order the one (i.e. tri-linear interpolation) instabilities are known to occur. To suppress those monotonicity constrains or flux-limiters need to be applies.
Exceptions
RuntimeErrorif inGrid is not aligned with mask or if its voxels are not uniform.
size_t getGrainSize ( ) const
inline
Returns
the grain-size used for multi-threading
Note
A grainsize of 0 implies serial execution
Scheme::SemiLagrangian getIntegrator ( ) const
inline

Return the integrator (see details in the table above)

Scheme::Limiter getLimiter ( ) const
inline

Retrun the limiter (see details above)

int getMaxDistance ( const VolumeGridT &  inGrid,
double  dt 
) const
inline
Returns
Returns the maximum distance in voxel units of inGrid that a particle can travel in the time-step dt when advected in the velocity field defined during construction.

This method is useful when dilating sparse volume grids to pad boundary regions. Excessive dilation can be computationally expensive so use this method to prevent or warn against run-away computation.

Exceptions
RuntimeErrorif inGrid does not have uniform voxels.
double getMaxVelocity ( ) const
inline

Return the maximum magnitude of the velocity in the advection velocity field defined during construction.

int getSubSteps ( ) const
inline
Returns
the number of sub-steps per integration (always larger than or equal to 1).
bool isLimiterOn ( ) const
inline

Return true if a limiter will be applied based on the current settings.

void setGrainSize ( size_t  grainsize)
inline

Set the grain-size used for multi-threading.

Note
A grainsize of 0 disables multi-threading
Warning
A small grainsize can degrade performance, both in terms of time and memory footprint!
void setIntegrator ( Scheme::SemiLagrangian  integrator)
inline

Set the integrator (see details in the table above)

void setLimiter ( Scheme::Limiter  limiter)
inline

Set the limiter (see details above)

void setSubSteps ( int  substeps)
inline

Set the number of sub-steps per integration.

Note
The only reason to increase the sub-step above its default value of one is to reduce the memory footprint due to significant dilation. Values smaller than 1 will be clamped to 1!
int spatialOrder ( ) const
inline

Return the spatial order of accuracy of the advection scheme.

Note
This is the optimal order in smooth regions. In non-smooth regions the flux-limiter will drop the order of accuracy to add numerical dissipation.
int temporalOrder ( ) const
inline

Return the temporal order of accuracy of the advection scheme.

Note
This is the optimal order in smooth regions. In non-smooth regions the flux-limiter will drop the order of accuracy to add numerical dissipation.