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...
|  | 
|  | VolumeAdvection (const VelocityGridT &velGrid, InterrupterType *interrupter=nullptr) | 
|  | Constructor. 
 | 
|  | 
| virtual | ~VolumeAdvection () | 
|  | 
| int | spatialOrder () const | 
|  | Return the spatial order of accuracy of the advection scheme. 
 | 
|  | 
| int | temporalOrder () const | 
|  | Return the temporal order of accuracy of the advection scheme. 
 | 
|  | 
| void | setIntegrator (Scheme::SemiLagrangian integrator) | 
|  | Set the integrator (see details in the table above) 
 | 
|  | 
| Scheme::SemiLagrangian | getIntegrator () const | 
|  | Return the integrator (see details in the table above) 
 | 
|  | 
| void | setLimiter (Scheme::Limiter limiter) | 
|  | Set the limiter (see details above) 
 | 
|  | 
| Scheme::Limiter | getLimiter () const | 
|  | Retrun the limiter (see details above) 
 | 
|  | 
| bool | isLimiterOn () const | 
|  | Return trueif a limiter will be applied based on the current settings.
 | 
|  | 
| size_t | getGrainSize () const | 
|  | 
| void | setGrainSize (size_t grainsize) | 
|  | Set the grain-size used for multi-threading. 
 | 
|  | 
| int | getSubSteps () const | 
|  | 
| void | setSubSteps (int substeps) | 
|  | Set the number of sub-steps per integration. 
 | 
|  | 
| double | getMaxVelocity () const | 
|  | Return the maximum magnitude of the velocity in the advection velocity field defined during construction. 
 | 
|  | 
| 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) | 
|  | 
template<typename VelocityGridT = Vec3fGrid, 
bool StaggeredVelocity = false, typename InterrupterType = util::NullInterrupter>
class openvdb::v12_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       |
///    ================================================================
///