HDK
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
VolumeAdvect.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (c) 2012-2017 DreamWorks Animation LLC
4 //
5 // All rights reserved. This software is distributed under the
6 // Mozilla Public License 2.0 ( http://www.mozilla.org/MPL/2.0/ )
7 //
8 // Redistributions of source code must retain the above copyright
9 // and license notice and the following restrictions and disclaimer.
10 //
11 // * Neither the name of DreamWorks Animation nor the names of
12 // its contributors may be used to endorse or promote products derived
13 // from this software without specific prior written permission.
14 //
15 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
16 // "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
17 // LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
18 // A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
19 // OWNER OR CONTRIBUTORS BE LIABLE FOR ANY INDIRECT, INCIDENTAL,
20 // SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
21 // LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
22 // DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
23 // THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
24 // (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
25 // OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
26 // IN NO EVENT SHALL THE COPYRIGHT HOLDERS' AND CONTRIBUTORS' AGGREGATE
27 // LIABILITY FOR ALL CLAIMS REGARDLESS OF THEIR BASIS EXCEED US$250.00.
28 //
29 ///////////////////////////////////////////////////////////////////////////
30 //
31 ///////////////////////////////////////////////////////////////////////////
32 //
33 /// @author Ken Museth
34 ///
35 /// @file VolumeAdvect.h
36 ///
37 /// @brief Sparse hyperbolic advection of volumes, e.g. a density or
38 /// velocity (vs a level set interface).
39 
40 #ifndef OPENVDB_TOOLS_VOLUME_ADVECT_HAS_BEEN_INCLUDED
41 #define OPENVDB_TOOLS_VOLUME_ADVECT_HAS_BEEN_INCLUDED
42 
43 #include <tbb/parallel_for.h>
44 #include <hboost/bind.hpp>
45 #include <hboost/function.hpp>
46 #include <openvdb/Types.h>
47 #include <openvdb/math/Math.h>
49 #include "Interpolation.h"// for Sampler
50 #include "VelocityFields.h" // for VelocityIntegrator
51 #include "Morphology.h"//for dilateActiveValues and dilateVoxels
52 #include "Prune.h"// for prune
53 #include "Statistics.h" // for extrema
54 
55 namespace openvdb {
57 namespace OPENVDB_VERSION_NAME {
58 namespace tools {
59 
60 
61 namespace Scheme {
62  /// @brief Numerical advections schemes.
64  /// @brief Flux-limiters employed to stabalize the second-order
65  /// advection schemes MacCormack and BFECC.
67 }
68 
69 /// @brief Performs advections of an arbitrary type of volume in a
70 /// static velocity field. The advections are performed by means
71 /// of various derivatives of Semi-Lagrangian integration, i.e.
72 /// backwards tracking along the hyperbolic characteristics
73 /// followed by interpolation.
74 ///
75 /// @note Optionally a limiter can be combined with the higher-order
76 /// integration schemes MacCormack and BFECC. There are two
77 /// types of limiters (CLAMP and REVERT) that supress
78 /// non-physical oscillations by means of either claminging or
79 /// reverting to a first-order schemes when the function is not
80 /// bounded by the cell values used for tri-linear interpolation.
81 ///
82 /// @verbatim The supported integrations schemes:
83 ///
84 /// ================================================================
85 /// | Lable | Accuracy | Integration Scheme | Interpolations |
86 /// | |Time/Space| | velocity/volume |
87 /// ================================================================
88 /// | SEMI | 1/1 | Semi-Lagrangian | 1/1 |
89 /// | MID | 2/1 | Mid-Point | 2/1 |
90 /// | RK3 | 3/1 | 3rd Order Runge-Kutta | 3/1 |
91 /// | RK4 | 4/1 | 4th Order Runge-Kutta | 4/1 |
92 /// | MAC | 2/2 | MacCormack | 2/2 |
93 /// | BFECC | 2/2 | BFECC | 3/2 |
94 /// ================================================================
95 /// @endverbatim
96 
97 template<typename VelocityGridT = Vec3fGrid,
98  bool StaggeredVelocity = false,
99  typename InterrupterType = util::NullInterrupter>
101 {
102 public:
103 
104  /// @brief Constructor
105  ///
106  /// @param velGrid Velocity grid responsible for the (passive) advection.
107  /// @param interrupter Optional interrupter used to prematurely end computations.
108  ///
109  /// @note The velocity field is assumed to be constant for the duration of the
110  /// advection.
111  VolumeAdvection(const VelocityGridT& velGrid, InterrupterType* interrupter = NULL)
112  : mVelGrid(velGrid)
113  , mInterrupter(interrupter)
114  , mIntegrator( Scheme::SEMI )
115  , mLimiter( Scheme::CLAMP )
116  , mGrainSize( 128 )
117  , mSubSteps( 1 )
118  {
119  math::Extrema e = extrema(velGrid.cbeginValueAll(), /*threading*/true);
120  e.add(velGrid.background().length());
121  mMaxVelocity = e.max();
122  }
123 
125  {
126  }
127 
128  /// @brief Return the spatial order of accuracy of the advection scheme
129  ///
130  /// @note This is the optimal order in smooth regions. In
131  /// non-smooth regions the flux-limiter will drop the order of
132  /// accuracy to add numerical dissipation.
133  int spatialOrder() const { return (mIntegrator == Scheme::MAC ||
134  mIntegrator == Scheme::BFECC) ? 2 : 1; }
135 
136  /// @brief Return the temporal order of accuracy of the advection scheme
137  ///
138  /// @note This is the optimal order in smooth regions. In
139  /// non-smooth regions the flux-limiter will drop the order of
140  /// accuracy to add numerical dissipation.
141  int temporalOrder() const {
142  switch (mIntegrator) {
143  case Scheme::SEMI: return 1;
144  case Scheme::MID: return 2;
145  case Scheme::RK3: return 3;
146  case Scheme::RK4: return 4;
147  case Scheme::BFECC:return 2;
148  case Scheme::MAC: return 2;
149  }
150  return 0;//should never reach this point
151  }
152 
153  /// @brief Set the integrator (see details in the table above)
154  void setIntegrator(Scheme::SemiLagrangian integrator) { mIntegrator = integrator; }
155 
156  /// @brief Return the integrator (see details in the table above)
157  Scheme::SemiLagrangian getIntegrator() const { return mIntegrator; }
158 
159  /// @brief Set the limiter (see details above)
160  void setLimiter(Scheme::Limiter limiter) { mLimiter = limiter; }
161 
162  /// @brief Retrun the limiter (see details above)
163  Scheme::Limiter getLimiter() const { return mLimiter; }
164 
165  /// @brief Return @c true if a limiter will be applied based on
166  /// the current settings.
167  bool isLimiterOn() const { return this->spatialOrder()>1 &&
168  mLimiter != Scheme::NO_LIMITER; }
169 
170  /// @return the grain-size used for multi-threading
171  /// @note A grainsize of 0 implies serial execution
172  size_t getGrainSize() const { return mGrainSize; }
173 
174  /// @brief Set the grain-size used for multi-threading
175  /// @note A grainsize of 0 disables multi-threading
176  /// @warning A small grainsize can degrade performance,
177  /// both in terms of time and memory footprint!
178  void setGrainSize(size_t grainsize) { mGrainSize = grainsize; }
179 
180  /// @return the number of sub-steps per integration (always larger
181  /// than or equal to 1).
182  int getSubSteps() const { return mSubSteps; }
183 
184  /// @brief Set the number of sub-steps per integration.
185  /// @note The only reason to increase the sub-step above its
186  /// default value of one is to reduce the memory footprint
187  /// due to significant dilation. Values smaller than 1 will
188  /// be clamped to 1!
189  void setSubSteps(int substeps) { mSubSteps = math::Max(1, substeps); }
190 
191  /// @brief Return the maximum magnitude of the velocity in the
192  /// advection velocity field defined during construction.
193  double getMaxVelocity() const { return mMaxVelocity; }
194 
195  /// @return Returns the maximum distance in voxel units of @a inGrid
196  /// that a particle can travel in the time-step @a dt when advected
197  /// in the velocity field defined during construction.
198  ///
199  /// @details This method is useful when dilating sparse volume
200  /// grids to pad boundary regions. Excessive dilation can be
201  /// computationally expensive so use this method to prevent
202  /// or warn against run-away computation.
203  ///
204  /// @throw RuntimeError if @a inGrid does not have uniform voxels.
205  template<typename VolumeGridT>
206  int getMaxDistance(const VolumeGridT& inGrid, double dt) const
207  {
208  if (!inGrid.hasUniformVoxels()) {
209  OPENVDB_THROW(RuntimeError, "Volume grid does not have uniform voxels!");
210  }
211  const double d = mMaxVelocity*math::Abs(dt)/inGrid.voxelSize()[0];
212  return static_cast<int>( math::RoundUp(d) );
213  }
214 
215  /// @return Returns a new grid that is the result of passive advection
216  /// of all the active values the input grid by @a timeStep.
217  ///
218  /// @param inGrid The input grid to be advected (unmodified)
219  /// @param timeStep Time-step of the Runge-Kutta integrator.
220  ///
221  /// @details This method will advect all of the active values in
222  /// the input @a inGrid. To achieve this a
223  /// deep-copy is dilated to account for the material
224  /// transport. This dilation step can be slow for large
225  /// time steps @a dt or a velocity field with large magnitudes.
226  ///
227  /// @warning If the VolumeSamplerT is of higher order than one
228  /// (i.e. tri-linear interpolation) instabilities are
229  /// known to occure. To suppress those monotonicity
230  /// constrains or flux-limiters need to be applies.
231  ///
232  /// @throw RuntimeError if @a inGrid does not have uniform voxels.
233  template<typename VolumeGridT,
234  typename VolumeSamplerT>//only C++11 allows for a default argument
235  typename VolumeGridT::Ptr advect(const VolumeGridT& inGrid, double timeStep)
236  {
237  typename VolumeGridT::Ptr outGrid = inGrid.deepCopy();
238  const double dt = timeStep/mSubSteps;
239  const int n = this->getMaxDistance(inGrid, dt);
240  dilateActiveValues( outGrid->tree(), n, NN_FACE, EXPAND_TILES);
241  this->template cook<VolumeGridT, VolumeSamplerT>(*outGrid, inGrid, dt);
242  for (int step = 1; step < mSubSteps; ++step) {
243  typename VolumeGridT::Ptr tmpGrid = outGrid->deepCopy();
244  dilateActiveValues( tmpGrid->tree(), n, NN_FACE, EXPAND_TILES);
245  this->template cook<VolumeGridT, VolumeSamplerT>(*tmpGrid, *outGrid, dt);
246  outGrid.swap( tmpGrid );
247  }
248 
249  return outGrid;
250  }
251 
252  /// @return Returns a new grid that is the result of
253  /// passive advection of the active values in @a inGrid
254  /// that intersect the active values in @c mask. The time
255  /// of the output grid is incremented by @a timeStep.
256  ///
257  /// @param inGrid The input grid to be advected (unmodified).
258  /// @param mask The mask of active values defining the active voxels
259  /// in @c inGrid on which to perform advection. Only
260  /// if a value is active in both grids will it be modified.
261  /// @param timeStep Time-step for a single Runge-Kutta integration step.
262  ///
263  ///
264  /// @details This method will advect all of the active values in
265  /// the input @a inGrid that intersects with the
266  /// active values in @a mask. To achieve this a
267  /// deep-copy is dilated to account for the material
268  /// transport and finally cropped to the intersection
269  /// with @a mask. The dilation step can be slow for large
270  /// time steps @a dt or fast moving velocity fields.
271  ///
272  /// @warning If the VolumeSamplerT is of higher order the one
273  /// (i.e. tri-linear interpolation) instabilities are
274  /// known to occure. To suppress those monotonicity
275  /// constrains or flux-limiters need to be applies.
276  ///
277  /// @throw RuntimeError if @a inGrid is not aligned with @a mask
278  /// or if its voxels are not uniform.
279  template<typename VolumeGridT,
280  typename MaskGridT,
281  typename VolumeSamplerT>//only C++11 allows for a default argument
282  typename VolumeGridT::Ptr advect(const VolumeGridT& inGrid, const MaskGridT& mask, double timeStep)
283  {
284  if (inGrid.transform() != mask.transform()) {
285  OPENVDB_THROW(RuntimeError, "Volume grid and mask grid are misaligned! Consider "
286  "resampling either of the two grids into the index space of the other.");
287  }
288  typename VolumeGridT::Ptr outGrid = inGrid.deepCopy();
289  const double dt = timeStep/mSubSteps;
290  const int n = this->getMaxDistance(inGrid, dt);
291  dilateActiveValues( outGrid->tree(), n, NN_FACE, EXPAND_TILES);
292  outGrid->topologyIntersection( mask );
293  pruneInactive( outGrid->tree(), mGrainSize>0, mGrainSize );
294  this->template cook<VolumeGridT, VolumeSamplerT>(*outGrid, inGrid, dt);
295  outGrid->topologyUnion( inGrid );
296 
297  for (int step = 1; step < mSubSteps; ++step) {
298  typename VolumeGridT::Ptr tmpGrid = outGrid->deepCopy();
299  dilateActiveValues( tmpGrid->tree(), n, NN_FACE, EXPAND_TILES);
300  tmpGrid->topologyIntersection( mask );
301  pruneInactive( tmpGrid->tree(), mGrainSize>0, mGrainSize );
302  this->template cook<VolumeGridT, VolumeSamplerT>(*tmpGrid, *outGrid, dt);
303  tmpGrid->topologyUnion( inGrid );
304  outGrid.swap( tmpGrid );
305  }
306  return outGrid;
307  }
308 
309 private:
310  // disallow copy construction and copy by assignment!
311  VolumeAdvection(const VolumeAdvection&);// not implemented
312  VolumeAdvection& operator=(const VolumeAdvection&);// not implemented
313 
314  void start(const char* str) const
315  {
316  if (mInterrupter) mInterrupter->start(str);
317  }
318  void stop() const
319  {
320  if (mInterrupter) mInterrupter->end();
321  }
322  bool interrupt() const
323  {
324  if (mInterrupter && util::wasInterrupted(mInterrupter)) {
325  tbb::task::self().cancel_group_execution();
326  return true;
327  }
328  return false;
329  }
330 
331  template<typename VolumeGridT, typename VolumeSamplerT>
332  void cook(VolumeGridT& outGrid, const VolumeGridT& inGrid, double dt)
333  {
334  switch (mIntegrator) {
335  case Scheme::SEMI: {
336  Advect<VolumeGridT, 1, VolumeSamplerT> adv(inGrid, *this);
337  adv.cook(outGrid, dt);
338  break;
339  }
340  case Scheme::MID: {
341  Advect<VolumeGridT, 2, VolumeSamplerT> adv(inGrid, *this);
342  adv.cook(outGrid, dt);
343  break;
344  }
345  case Scheme::RK3: {
346  Advect<VolumeGridT, 3, VolumeSamplerT> adv(inGrid, *this);
347  adv.cook(outGrid, dt);
348  break;
349  }
350  case Scheme::RK4: {
351  Advect<VolumeGridT, 4, VolumeSamplerT> adv(inGrid, *this);
352  adv.cook(outGrid, dt);
353  break;
354  }
355  case Scheme::BFECC: {
356  Advect<VolumeGridT, 1, VolumeSamplerT> adv(inGrid, *this);
357  adv.cook(outGrid, dt);
358  break;
359  }
360  case Scheme::MAC: {
361  Advect<VolumeGridT, 1, VolumeSamplerT> adv(inGrid, *this);
362  adv.cook(outGrid, dt);
363  break;
364  }
365  default:
366  OPENVDB_THROW(ValueError, "Spatial difference scheme not supported!");
367  }
368  pruneInactive(outGrid.tree(), mGrainSize>0, mGrainSize);
369  }
370 
371  // Private class that implements the multi-threaded advection
372  template<typename VolumeGridT, size_t OrderRK, typename SamplerT> struct Advect;
373 
374  // Private member data of VolumeAdvection
375  const VelocityGridT& mVelGrid;
376  double mMaxVelocity;
377  InterrupterType* mInterrupter;
378  Scheme::SemiLagrangian mIntegrator;
379  Scheme::Limiter mLimiter;
380  size_t mGrainSize;
381  int mSubSteps;
382 };//end of VolumeAdvection class
383 
384 // Private class that implements the multi-threaded advection
385 template<typename VelocityGridT, bool StaggeredVelocity, typename InterrupterType>
386 template<typename VolumeGridT, size_t OrderRK, typename SamplerT>
387 struct VolumeAdvection<VelocityGridT, StaggeredVelocity, InterrupterType>::Advect
388 {
389  typedef typename VolumeGridT::TreeType TreeT;
390  typedef typename VolumeGridT::ConstAccessor AccT;
391  typedef typename TreeT::ValueType ValueT;
392  typedef typename tree::LeafManager<TreeT> LeafManagerT;
393  typedef typename LeafManagerT::LeafNodeType LeafNodeT;
394  typedef typename LeafManagerT::LeafRange LeafRangeT;
395  typedef VelocityIntegrator<VelocityGridT, StaggeredVelocity> VelocityIntegratorT;
396  typedef typename VelocityIntegratorT::ElementType RealT;
397  typedef typename TreeT::LeafNodeType::ValueOnIter VoxelIterT;
398 
399  Advect(const VolumeGridT& inGrid, const VolumeAdvection& parent)
400  : mTask(0)
401  , mInGrid(&inGrid)
402  , mVelocityInt(parent.mVelGrid)
403  , mParent(&parent)
404  {
405  }
406  inline void cook(const LeafRangeT& range)
407  {
408  if (mParent->mGrainSize > 0) {
409  tbb::parallel_for(range, *this);
410  } else {
411  (*this)(range);
412  }
413  }
414  void operator()(const LeafRangeT& range) const
415  {
416  assert(mTask);
417  mTask(const_cast<Advect*>(this), range);
418  }
419  void cook(VolumeGridT& outGrid, double time_step)
420  {
421  mParent->start("Advecting volume");
422  LeafManagerT manager(outGrid.tree(), mParent->spatialOrder()==2 ? 1 : 0);
423  const LeafRangeT range = manager.leafRange(mParent->mGrainSize);
424  const RealT dt = static_cast<RealT>(-time_step);//method of characteristics backtracks
425  if (mParent->mIntegrator == Scheme::MAC) {
426  mTask = hboost::bind(&Advect::rk, hboost::placeholders::_1, hboost::placeholders::_2, dt, 0, mInGrid);//out[0]=forward
427  this->cook(range);
428  mTask = hboost::bind(&Advect::rk, hboost::placeholders::_1, hboost::placeholders::_2,-dt, 1, &outGrid);//out[1]=backward
429  this->cook(range);
430  mTask = hboost::bind(&Advect::mac, hboost::placeholders::_1, hboost::placeholders::_2);//out[0] = out[0] + (in[0] - out[1])/2
431  this->cook(range);
432  } else if (mParent->mIntegrator == Scheme::BFECC) {
433  mTask = hboost::bind(&Advect::rk, hboost::placeholders::_1, hboost::placeholders::_2, dt, 0, mInGrid);//out[0]=forward
434  this->cook(range);
435  mTask = hboost::bind(&Advect::rk, hboost::placeholders::_1, hboost::placeholders::_2,-dt, 1, &outGrid);//out[1]=backward
436  this->cook(range);
437  mTask = hboost::bind(&Advect::bfecc, hboost::placeholders::_1, hboost::placeholders::_2);//out[0] = (3*in[0] - out[1])/2
438  this->cook(range);
439  mTask = hboost::bind(&Advect::rk, hboost::placeholders::_1, hboost::placeholders::_2, dt, 1, &outGrid);//out[1]=forward
440  this->cook(range);
441  manager.swapLeafBuffer(1);// out[0] = out[1]
442  } else {// SEMI, MID, RK3 and RK4
443  mTask = hboost::bind(&Advect::rk, hboost::placeholders::_1, hboost::placeholders::_2, dt, 0, mInGrid);//forward
444  this->cook(range);
445  }
446 
447  if (mParent->spatialOrder()==2) manager.removeAuxBuffers();
448 
449  mTask = hboost::bind(&Advect::limiter, hboost::placeholders::_1, hboost::placeholders::_2, dt);// out[0] = limiter( out[0] )
450  this->cook(range);
451 
452  mParent->stop();
453  }
454  // Last step of the MacCormack scheme: out[0] = out[0] + (in[0] - out[1])/2
455  void mac(const LeafRangeT& range) const
456  {
457  if (mParent->interrupt()) return;
458  assert( mParent->mIntegrator == Scheme::MAC );
459  AccT acc = mInGrid->getAccessor();
460  for (typename LeafRangeT::Iterator leafIter = range.begin(); leafIter; ++leafIter) {
461  ValueT* out0 = leafIter.buffer( 0 ).data();// forward
462  const ValueT* out1 = leafIter.buffer( 1 ).data();// backward
463  const LeafNodeT* leaf = acc.probeConstLeaf( leafIter->origin() );
464  if (leaf !=NULL) {
465  const ValueT* in0 = leaf->buffer().data();
466  for (VoxelIterT voxelIter = leafIter->beginValueOn(); voxelIter; ++voxelIter) {
467  const Index i = voxelIter.pos();
468  out0[i] += RealT(0.5) * ( in0[i] - out1[i] );
469  }
470  } else {
471  for (VoxelIterT voxelIter = leafIter->beginValueOn(); voxelIter; ++voxelIter) {
472  const Index i = voxelIter.pos();
473  out0[i] += RealT(0.5) * ( acc.getValue(voxelIter.getCoord()) - out1[i] );
474  }//loop over active voxels
475  }
476  }//loop over leaf nodes
477  }
478  // Intermediate step in the BFECC scheme: out[0] = (3*in[0] - out[1])/2
479  void bfecc(const LeafRangeT& range) const
480  {
481  if (mParent->interrupt()) return;
482  assert( mParent->mIntegrator == Scheme::BFECC );
483  AccT acc = mInGrid->getAccessor();
484  for (typename LeafRangeT::Iterator leafIter = range.begin(); leafIter; ++leafIter) {
485  ValueT* out0 = leafIter.buffer( 0 ).data();// forward
486  const ValueT* out1 = leafIter.buffer( 1 ).data();// backward
487  const LeafNodeT* leaf = acc.probeConstLeaf(leafIter->origin());
488  if (leaf !=NULL) {
489  const ValueT* in0 = leaf->buffer().data();
490  for (VoxelIterT voxelIter = leafIter->beginValueOn(); voxelIter; ++voxelIter) {
491  const Index i = voxelIter.pos();
492  out0[i] = RealT(0.5)*( RealT(3)*in0[i] - out1[i] );
493  }//loop over active voxels
494  } else {
495  for (VoxelIterT voxelIter = leafIter->beginValueOn(); voxelIter; ++voxelIter) {
496  const Index i = voxelIter.pos();
497  out0[i] = RealT(0.5)*( RealT(3)*acc.getValue(voxelIter.getCoord()) - out1[i] );
498  }//loop over active voxels
499  }
500  }//loop over leaf nodes
501  }
502  // Semi-Lagrangian integration with Runge-Kutta of various orders (1->4)
503  void rk(const LeafRangeT& range, RealT dt, size_t n, const VolumeGridT* grid) const
504  {
505  if (mParent->interrupt()) return;
506  const math::Transform& xform = mInGrid->transform();
507  AccT acc = grid->getAccessor();
508  for (typename LeafRangeT::Iterator leafIter = range.begin(); leafIter; ++leafIter) {
509  ValueT* phi = leafIter.buffer( n ).data();
510  for (VoxelIterT voxelIter = leafIter->beginValueOn(); voxelIter; ++voxelIter) {
511  ValueT& value = phi[voxelIter.pos()];
512  Vec3d wPos = xform.indexToWorld(voxelIter.getCoord());
513  mVelocityInt.template rungeKutta<OrderRK, Vec3d>(dt, wPos);
514  value = SamplerT::sample(acc, xform.worldToIndex(wPos));
515  }//loop over active voxels
516  }//loop over leaf nodes
517  }
518  void limiter(const LeafRangeT& range, RealT dt) const
519  {
520  if (mParent->interrupt()) return;
521  const bool doLimiter = mParent->isLimiterOn();
522  const bool doClamp = mParent->mLimiter == Scheme::CLAMP;
523  ValueT data[2][2][2], vMin, vMax;
524  const math::Transform& xform = mInGrid->transform();
525  AccT acc = mInGrid->getAccessor();
526  const ValueT backg = mInGrid->background();
527  for (typename LeafRangeT::Iterator leafIter = range.begin(); leafIter; ++leafIter) {
528  ValueT* phi = leafIter.buffer( 0 ).data();
529  for (VoxelIterT voxelIter = leafIter->beginValueOn(); voxelIter; ++voxelIter) {
530  ValueT& value = phi[voxelIter.pos()];
531 
532  if ( doLimiter ) {
533  assert(OrderRK == 1);
534  Vec3d wPos = xform.indexToWorld(voxelIter.getCoord());
535  mVelocityInt.template rungeKutta<1, Vec3d>(dt, wPos);// Explicit Euler
536  Vec3d iPos = xform.worldToIndex(wPos);
537  Coord ijk = Coord::floor( iPos );
538  BoxSampler::getValues(data, acc, ijk);
539  BoxSampler::extrema(data, vMin, vMax);
540  if ( doClamp ) {
541  value = math::Clamp( value, vMin, vMax);
542  } else if (value < vMin || value > vMax ) {
543  iPos -= Vec3R(ijk[0], ijk[1], ijk[2]);//unit coordinates
544  value = BoxSampler::trilinearInterpolation( data, iPos );
545  }
546  }
547 
548  if (math::isApproxEqual(value, backg, math::Delta<ValueT>::value())) {
549  value = backg;
550  leafIter->setValueOff( voxelIter.pos() );
551  }
552  }//loop over active voxels
553  }//loop over leaf nodes
554  }
555  // Public member data of the private Advect class
556 
557  typename hboost::function<void (Advect*, const LeafRangeT&)> mTask;
558  const VolumeGridT* mInGrid;
559  const VelocityIntegratorT mVelocityInt;// lightweight!
560  const VolumeAdvection* mParent;
561 };// end of private member class Advect
562 
563 } // namespace tools
564 } // namespace OPENVDB_VERSION_NAME
565 } // namespace openvdb
566 
567 #endif // OPENVDB_TOOLS_VOLUME_ADVECT_HAS_BEEN_INCLUDED
568 
569 // Copyright (c) 2012-2017 DreamWorks Animation LLC
570 // All rights reserved. This software is distributed under the
571 // Mozilla Public License 2.0 ( http://www.mozilla.org/MPL/2.0/ )
GLenum GLint * range
Definition: glcorearb.h:1924
math::Vec3< Real > Vec3R
Definition: Types.h:75
SemiLagrangian
Numerical advections schemes.
Definition: VolumeAdvect.h:63
int getMaxDistance(const VolumeGridT &inGrid, double dt) const
Definition: VolumeAdvect.h:206
void setGrainSize(size_t grainsize)
Set the grain-size used for multi-threading.
Definition: VolumeAdvect.h:178
GLuint start
Definition: glcorearb.h:474
static ValueT trilinearInterpolation(ValueT(&data)[N][N][N], const Vec3R &uvw)
int temporalOrder() const
Return the temporal order of accuracy of the advection scheme.
Definition: VolumeAdvect.h:141
void setSubSteps(int substeps)
Set the number of sub-steps per integration.
Definition: VolumeAdvect.h:189
GLint GLuint mask
Definition: glcorearb.h:123
static void extrema(ValueT(&data)[N][N][N], ValueT &vMin, ValueT &vMax)
Find the minimum and maximum values of the eight cell values in @ data.
Dummy NOOP interrupter class defining interface.
png_uint_32 i
Definition: png.h:2877
float RoundUp(float x)
Return x rounded up to the nearest integer.
Definition: Math.h:753
Performs advections of an arbitrary type of volume in a static velocity field. The advections are per...
Definition: VolumeAdvect.h:100
VolumeGridT::Ptr advect(const VolumeGridT &inGrid, const MaskGridT &mask, double timeStep)
Definition: VolumeAdvect.h:282
Limiter
Flux-limiters employed to stabalize the second-order advection schemes MacCormack and BFECC...
Definition: VolumeAdvect.h:66
GLdouble n
Definition: glcorearb.h:2007
Coord Abs(const Coord &xyz)
Definition: Coord.h:254
Scheme::SemiLagrangian getIntegrator() const
Return the integrator (see details in the table above)
Definition: VolumeAdvect.h:157
#define OPENVDB_VERSION_NAME
Definition: version.h:43
General-purpose arithmetic and comparison routines, most of which accept arbitrary value types (or at...
VolumeAdvection(const VelocityGridT &velGrid, InterrupterType *interrupter=NULL)
Constructor.
Definition: VolumeAdvect.h:111
Defined various multi-threaded utility functions for trees.
GLboolean * data
Definition: glcorearb.h:130
static void getValues(ValueT(&data)[N][N][N], const TreeT &inTree, Coord ijk)
Import all eight values from inTree to support tri-linear interpolation.
bool isApproxEqual(const Type &a, const Type &b)
Return true if a is equal to b to within the default floating-point comparison tolerance.
Definition: Math.h:370
int floor(T x)
Definition: ImathFun.h:150
Type Clamp(Type x, Type min, Type max)
Return x clamped to [min, max].
Definition: Math.h:246
bool isLimiterOn() const
Return true if a limiter will be applied based on the current settings.
Definition: VolumeAdvect.h:167
double getMaxVelocity() const
Return the maximum magnitude of the velocity in the advection velocity field defined during construct...
Definition: VolumeAdvect.h:193
GLsizei const GLfloat * value
Definition: glcorearb.h:823
void setIntegrator(Scheme::SemiLagrangian integrator)
Set the integrator (see details in the table above)
Definition: VolumeAdvect.h:154
Defines two simple wrapper classes for advection velocity fields as well as VelocitySampler and Veloc...
int spatialOrder() const
Return the spatial order of accuracy of the advection scheme.
Definition: VolumeAdvect.h:133
VolumeGridT::Ptr advect(const VolumeGridT &inGrid, double timeStep)
Definition: VolumeAdvect.h:235
Implementation of morphological dilation and erosion.
Functions to efficiently compute histograms, extremas (min/max) and statistics (mean, variance, etc.) of grid values.
Scheme::Limiter getLimiter() const
Retrun the limiter (see details above)
Definition: VolumeAdvect.h:163
#define OPENVDB_USE_VERSION_NAMESPACE
Definition: version.h:71
bool wasInterrupted(T *i, int percent=-1)
void setLimiter(Scheme::Limiter limiter)
Set the limiter (see details above)
Definition: VolumeAdvect.h:160
const Type & Max(const Type &a, const Type &b)
Return the maximum of two values.
Definition: Math.h:561
double max() const
Return the maximum value.
Definition: Stats.h:141
OPENVDB_STATIC_SPECIALIZATION void dilateActiveValues(TreeType &tree, int iterations=1, NearestNeighbors nn=NN_FACE, TilePolicy mode=PRESERVE_TILES)
Topologically dilate all active values (i.e. both voxels and tiles) in a tree using one of three near...
Definition: Morphology.h:1074
math::Extrema extrema(const IterT &iter, bool threaded=true)
Iterate over a scalar grid and compute extrema (min/max) of the values of the voxels that are visited...
Definition: Statistics.h:379
void add(double val)
Add a single sample.
Definition: Stats.h:119
void pruneInactive(TreeT &tree, bool threaded=true, size_t grainSize=1)
Reduce the memory footprint of a tree by replacing with background tiles any nodes whose values are a...
Definition: Prune.h:394
#define OPENVDB_THROW(exception, message)
Definition: Exceptions.h:101
This class computes the minimum and maximum values of a population of floating-point values...
Definition: Stats.h:105