HDK
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ParticlesToLevelSet.h
Go to the documentation of this file.
1 // Copyright Contributors to the OpenVDB Project
2 // SPDX-License-Identifier: MPL-2.0
3 
4 /// @author Ken Museth
5 ///
6 /// @file tools/ParticlesToLevelSet.h
7 ///
8 /// @brief Rasterize particles with position, radius and velocity
9 /// into either a boolean mask grid or a narrow-band level set grid.
10 ///
11 /// @details Optionally, arbitrary attributes on the particles can be transferred,
12 /// resulting in additional output grids with the same topology as the main grid.
13 ///
14 /// @note Particle to level set conversion is intended to be combined with
15 /// some kind of surface postprocessing, using
16 /// @vdblink::tools::LevelSetFilter LevelSetFilter@endlink, for example.
17 /// Without such postprocessing the generated surface is typically too noisy and blobby.
18 /// However, it serves as a great and fast starting point for subsequent
19 /// level set surface processing and convolution.
20 ///
21 /// @details For particle access, any class with the following interface may be used
22 /// (see the unit test or the From Particles Houdini SOP for practical examples):
23 /// @code
24 /// struct ParticleList
25 /// {
26 /// // Return the total number of particles in the list.
27 /// // Always required!
28 /// size_t size() const;
29 ///
30 /// // Get the world-space position of the nth particle.
31 /// // Required by rasterizeSpheres().
32 /// void getPos(size_t n, Vec3R& xyz) const;
33 ///
34 /// // Get the world-space position and radius of the nth particle.
35 /// // Required by rasterizeSpheres().
36 /// void getPosRad(size_t n, Vec3R& xyz, Real& radius) const;
37 ///
38 /// // Get the world-space position, radius and velocity of the nth particle.
39 /// // Required by rasterizeTrails().
40 /// void getPosRadVel(size_t n, Vec3R& xyz, Real& radius, Vec3R& velocity) const;
41 ///
42 /// // Get the value of the nth particle's user-defined attribute (of type @c AttributeType).
43 /// // Required only if attribute transfer is enabled in ParticlesToLevelSet.
44 /// void getAtt(size_t n, AttributeType& att) const;
45 /// };
46 /// @endcode
47 ///
48 /// Some functions accept an interrupter argument. This refers to any class
49 /// with the following interface:
50 /// @code
51 /// struct Interrupter
52 /// {
53 /// void start(const char* name = nullptr) // called when computations begin
54 /// void end() // called when computations end
55 /// bool wasInterrupted(int percent=-1) // return true to abort computation
56 /// };
57 /// @endcode
58 ///
59 /// The default interrupter is @vdblink::util::NullInterrupter NullInterrupter@endlink,
60 /// for which all calls are no-ops that incur no computational overhead.
61 
62 #ifndef OPENVDB_TOOLS_PARTICLES_TO_LEVELSET_HAS_BEEN_INCLUDED
63 #define OPENVDB_TOOLS_PARTICLES_TO_LEVELSET_HAS_BEEN_INCLUDED
64 
65 #include <openvdb/Types.h>
66 #include <openvdb/Grid.h>
67 #include <openvdb/math/Math.h>
68 #include <openvdb/math/Transform.h>
70 #include <openvdb/util/logging.h>
73 
74 #include "Composite.h" // for csgUnion()
75 #include "PointPartitioner.h"
76 #include "Prune.h"
77 #include "SignedFloodFill.h"
78 
79 #include <tbb/parallel_reduce.h>
80 #include <tbb/blocked_range.h>
81 
82 #include <functional>
83 #include <iostream>
84 #include <type_traits>
85 #include <vector>
86 
87 
88 namespace openvdb {
90 namespace OPENVDB_VERSION_NAME {
91 namespace tools {
92 
93 /// @brief Populate a scalar, floating-point grid with CSG-unioned level set spheres
94 /// described by the given particle positions and radii.
95 /// @details For more control over the output, including attribute transfer,
96 /// use the ParticlesToLevelSet class directly.
97 template<typename GridT, typename ParticleListT, typename InterrupterT = util::NullInterrupter>
98 inline void particlesToSdf(const ParticleListT&, GridT&, InterrupterT* = nullptr);
99 
100 /// @brief Populate a scalar, floating-point grid with fixed-size, CSG-unioned
101 /// level set spheres described by the given particle positions and the specified radius.
102 /// @details For more control over the output, including attribute transfer,
103 /// use the ParticlesToLevelSet class directly.
104 template<typename GridT, typename ParticleListT, typename InterrupterT = util::NullInterrupter>
105 inline void particlesToSdf(const ParticleListT&, GridT&, Real radius, InterrupterT* = nullptr);
106 
107 /// @brief Populate a scalar, floating-point grid with CSG-unioned trails
108 /// of level set spheres with decreasing radius, where the starting position and radius
109 /// and the direction of each trail is given by particle attributes.
110 /// @details For more control over the output, including attribute transfer,
111 /// use the ParticlesToLevelSet class directly.
112 /// @note The @a delta parameter controls the distance between spheres in a trail.
113 /// Be careful not to use too small a value.
114 template<typename GridT, typename ParticleListT, typename InterrupterT = util::NullInterrupter>
115 inline void particleTrailsToSdf(const ParticleListT&, GridT&, Real delta=1, InterrupterT* =nullptr);
116 
117 /// @brief Activate a boolean grid wherever it intersects the spheres
118 /// described by the given particle positions and radii.
119 /// @details For more control over the output, including attribute transfer,
120 /// use the ParticlesToLevelSet class directly.
121 template<typename GridT, typename ParticleListT, typename InterrupterT = util::NullInterrupter>
122 inline void particlesToMask(const ParticleListT&, GridT&, InterrupterT* = nullptr);
123 
124 /// @brief Activate a boolean grid wherever it intersects the fixed-size spheres
125 /// described by the given particle positions and the specified radius.
126 /// @details For more control over the output, including attribute transfer,
127 /// use the ParticlesToLevelSet class directly.
128 template<typename GridT, typename ParticleListT, typename InterrupterT = util::NullInterrupter>
129 inline void particlesToMask(const ParticleListT&, GridT&, Real radius, InterrupterT* = nullptr);
130 
131 /// @brief Activate a boolean grid wherever it intersects trails of spheres
132 /// with decreasing radius, where the starting position and radius and the direction
133 /// of each trail is given by particle attributes.
134 /// @details For more control over the output, including attribute transfer,
135 /// use the ParticlesToLevelSet class directly.
136 /// @note The @a delta parameter controls the distance between spheres in a trail.
137 /// Be careful not to use too small a value.
138 template<typename GridT, typename ParticleListT, typename InterrupterT = util::NullInterrupter>
139 inline void particleTrailsToMask(const ParticleListT&, GridT&,Real delta=1,InterrupterT* =nullptr);
140 
141 
142 ////////////////////////////////////////
143 
144 /// @cond OPENVDB_DOCS_INTERNAL
145 
146 namespace p2ls_internal {
147 // This is a simple type that combines a distance value and a particle
148 // attribute. It's required for attribute transfer which is performed
149 // in the ParticlesToLevelSet::Raster member class defined below.
150 /// @private
151 template<typename VisibleT, typename BlindT> class BlindData;
152 }
153 
154 /// @endcond
155 
156 template<typename SdfGridT,
157  typename AttributeT = void,
158  typename InterrupterT = util::NullInterrupter>
160 {
161 public:
163  using InterrupterType = InterrupterT;
164 
165  using SdfGridType = SdfGridT;
166  using SdfType = typename SdfGridT::ValueType;
167 
169  using AttGridType = typename SdfGridT::template ValueConverter<AttType>::Type;
170 
172 
173  /// @brief Constructor using an existing boolean or narrow-band level set grid
174  ///
175  /// @param grid grid into which particles are rasterized
176  /// @param interrupt callback to interrupt a long-running process
177  ///
178  /// @details If the input grid is already populated with signed distances,
179  /// particles are unioned onto the existing level set surface.
180  ///
181  /// @details The width in voxel units of the generated narrow band level set
182  /// is given by 2&times;<I>background</I>/<I>dx</I>, where @a background
183  /// is the background value stored in the grid and @a dx is the voxel size
184  /// derived from the transform associated with the grid.
185  /// Also note that &minus;<I>background</I> corresponds to the constant value
186  /// inside the generated narrow-band level set.
187  ///
188  /// @note If attribute transfer is enabled, i.e., if @c AttributeT is not @c void,
189  /// attributes are generated only for voxels that overlap with particles,
190  /// not for any other preexisting voxels (for which no attributes exist!).
191  explicit ParticlesToLevelSet(SdfGridT& grid, InterrupterT* interrupt = nullptr);
192 
193  ~ParticlesToLevelSet() { delete mBlindGrid; }
194 
195  /// @brief This method syncs up the level set and attribute grids
196  /// and therefore needs to be called before any of those grids are
197  /// used and after the last call to any of the rasterizer methods.
198  /// @details It has no effect or overhead if attribute transfer is disabled
199  /// (i.e., if @c AttributeT is @c void) and @a prune is @c false.
200  ///
201  /// @note Avoid calling this method more than once, and call it only after
202  /// all the particles have been rasterized.
203  void finalize(bool prune = false);
204 
205  /// @brief Return a pointer to the grid containing the optional user-defined attribute.
206  /// @warning If attribute transfer is disabled (i.e., if @c AttributeT is @c void)
207  /// or if @link finalize() finalize@endlink is not called, the pointer will be null.
208  typename AttGridType::Ptr attributeGrid() { return mAttGrid; }
209 
210  /// @brief Return the size of a voxel in world units.
211  Real getVoxelSize() const { return mDx; }
212 
213  /// @brief Return the half-width of the narrow band in voxel units.
214  Real getHalfWidth() const { return mHalfWidth; }
215 
216  /// @brief Return the smallest radius allowed in voxel units.
217  Real getRmin() const { return mRmin; }
218  /// @brief Set the smallest radius allowed in voxel units.
219  void setRmin(Real Rmin) { mRmin = math::Max(Real(0),Rmin); }
220 
221  /// @brief Return the largest radius allowed in voxel units.
222  Real getRmax() const { return mRmax; }
223  /// @brief Set the largest radius allowed in voxel units.
224  void setRmax(Real Rmax) { mRmax = math::Max(mRmin,Rmax); }
225 
226  /// @brief Return @c true if any particles were ignored due to their size.
227  bool ignoredParticles() const { return mMinCount>0 || mMaxCount>0; }
228  /// @brief Return the number of particles that were ignored because they were
229  /// smaller than the minimum radius.
230  size_t getMinCount() const { return mMinCount; }
231  /// @brief Return the number of particles that were ignored because they were
232  /// larger than the maximum radius.
233  size_t getMaxCount() const { return mMaxCount; }
234 
235  /// @brief Return the grain size used for threading
236  int getGrainSize() const { return mGrainSize; }
237  /// @brief Set the grain size used for threading.
238  /// @note A grain size of zero or less disables threading.
239  void setGrainSize(int grainSize) { mGrainSize = grainSize; }
240 
241  /// @brief Rasterize each particle as a sphere with the particle's position and radius.
242  /// @details For level set output, all spheres are CSG-unioned.
243  template<typename ParticleListT>
244  void rasterizeSpheres(const ParticleListT& pa);
245 
246  /// @brief Rasterize each particle as a sphere with the particle's position
247  /// and a fixed radius.
248  /// @details For level set output, all spheres are CSG-unioned.
249  ///
250  /// @param pa particles with positions
251  /// @param radius fixed sphere radius in world units.
252  template<typename ParticleListT>
253  void rasterizeSpheres(const ParticleListT& pa, Real radius);
254 
255  /// @brief Rasterize each particle as a trail comprising the CSG union
256  /// of spheres of decreasing radius.
257  ///
258  /// @param pa particles with position, radius and velocity.
259  /// @param delta controls the distance between sphere instances
260  ///
261  /// @warning Be careful not to use too small values for @a delta,
262  /// since this can lead to excessive computation per trail (which the
263  /// interrupter can't stop).
264  ///
265  /// @note The direction of a trail is opposite to that of the velocity vector,
266  /// and its length is given by the magnitude of the velocity.
267  /// The radius at the head of the trail is given by the radius of the particle,
268  /// and the radius at the tail is @a Rmin voxel units, which has
269  /// a default value of 1.5 corresponding to the Nyquist frequency!
270  template<typename ParticleListT>
271  void rasterizeTrails(const ParticleListT& pa, Real delta=1.0);
272 
273 private:
274  using BlindType = p2ls_internal::BlindData<SdfType, AttType>;
275  using BlindGridType = typename SdfGridT::template ValueConverter<BlindType>::Type;
276 
277  /// Class with multi-threaded implementation of particle rasterization
278  template<typename ParticleListT, typename GridT> struct Raster;
279 
280  SdfGridType* mSdfGrid;
281  typename AttGridType::Ptr mAttGrid;
282  BlindGridType* mBlindGrid;
283  InterrupterT* mInterrupter;
284  Real mDx, mHalfWidth;
285  Real mRmin, mRmax; // ignore particles outside this range of radii in voxel
286  size_t mMinCount, mMaxCount; // counters for ignored particles
287  int mGrainSize;
288 }; // class ParticlesToLevelSet
289 
290 
291 template<typename SdfGridT, typename AttributeT, typename InterrupterT>
293 ParticlesToLevelSet(SdfGridT& grid, InterrupterT* interrupter) :
294  mSdfGrid(&grid),
295  mBlindGrid(nullptr),
296  mInterrupter(interrupter),
297  mDx(grid.voxelSize()[0]),
298  mHalfWidth(grid.background()/mDx),
299  mRmin(1.5),// corresponds to the Nyquist grid sampling frequency
300  mRmax(100.0),// corresponds to a huge particle (probably too large!)
301  mMinCount(0),
302  mMaxCount(0),
303  mGrainSize(1)
304 {
305  if (!mSdfGrid->hasUniformVoxels()) {
306  OPENVDB_THROW(RuntimeError, "ParticlesToLevelSet only supports uniform voxels!");
307  }
308  if (!DisableT::value) {
309  mBlindGrid = new BlindGridType(BlindType(grid.background()));
310  mBlindGrid->setTransform(mSdfGrid->transform().copy());
311  }
312 }
313 
314 template<typename SdfGridT, typename AttributeT, typename InterrupterT>
315 template<typename ParticleListT>
317 rasterizeSpheres(const ParticleListT& pa)
318 {
319  if (DisableT::value) {
320  Raster<ParticleListT, SdfGridT> r(*this, mSdfGrid, pa);
321  r.rasterizeSpheres();
322  } else {
323  Raster<ParticleListT, BlindGridType> r(*this, mBlindGrid, pa);
324  r.rasterizeSpheres();
325  }
326 }
327 
328 template<typename SdfGridT, typename AttributeT, typename InterrupterT>
329 template<typename ParticleListT>
331 rasterizeSpheres(const ParticleListT& pa, Real radius)
332 {
333  if (DisableT::value) {
334  Raster<ParticleListT, SdfGridT> r(*this, mSdfGrid, pa);
335  r.rasterizeSpheres(radius/mDx);
336  } else {
337  Raster<ParticleListT, BlindGridType> r(*this, mBlindGrid, pa);
338  r.rasterizeSpheres(radius/mDx);
339  }
340 }
341 
342 template<typename SdfGridT, typename AttributeT, typename InterrupterT>
343 template<typename ParticleListT>
345 rasterizeTrails(const ParticleListT& pa, Real delta)
346 {
347  if (DisableT::value) {
348  Raster<ParticleListT, SdfGridT> r(*this, mSdfGrid, pa);
349  r.rasterizeTrails(delta);
350  } else {
351  Raster<ParticleListT, BlindGridType> r(*this, mBlindGrid, pa);
352  r.rasterizeTrails(delta);
353  }
354 }
355 
356 
357 template<typename SdfGridT, typename AttributeT, typename InterrupterT>
358 inline void
360 {
362 
363  if (!mBlindGrid) {
364  if (prune) {
365  if (OutputIsMask) {
366  tools::prune(mSdfGrid->tree());
367  } else {
368  tools::pruneLevelSet(mSdfGrid->tree());
369  }
370  }
371  return;
372  }
373 
374  if (prune) tools::prune(mBlindGrid->tree());
375 
376  using AttTreeT = typename AttGridType::TreeType;
377  using AttLeafT = typename AttTreeT::LeafNodeType;
378  using BlindTreeT = typename BlindGridType::TreeType;
379  using BlindLeafIterT = typename BlindTreeT::LeafCIter;
380  using BlindLeafT = typename BlindTreeT::LeafNodeType;
381  using SdfTreeT = typename SdfGridType::TreeType;
382  using SdfLeafT = typename SdfTreeT::LeafNodeType;
383 
384  // Use topology copy constructors since output grids have the same topology as mBlindDataGrid
385  const BlindTreeT& blindTree = mBlindGrid->tree();
386 
387  // Create the output attribute grid.
388  typename AttTreeT::Ptr attTree(new AttTreeT(
389  blindTree, blindTree.background().blind(), openvdb::TopologyCopy()));
390  // Note this overwrites any existing attribute grids!
391  mAttGrid = typename AttGridType::Ptr(new AttGridType(attTree));
392  mAttGrid->setTransform(mBlindGrid->transform().copy());
393 
394  typename SdfTreeT::Ptr sdfTree; // the output mask or level set tree
395 
396  // Extract the attribute grid and the mask or level set grid from mBlindDataGrid.
397  if (OutputIsMask) {
398  sdfTree.reset(new SdfTreeT(blindTree,
399  /*off=*/SdfType(0), /*on=*/SdfType(1), TopologyCopy()));
400 
401  // Copy leaf voxels in parallel.
402  tree::LeafManager<AttTreeT> leafNodes(*attTree);
403  leafNodes.foreach([&](AttLeafT& attLeaf, size_t /*leafIndex*/) {
404  if (const auto* blindLeaf = blindTree.probeConstLeaf(attLeaf.origin())) {
405  for (auto iter = attLeaf.beginValueOn(); iter; ++iter) {
406  const auto pos = iter.pos();
407  attLeaf.setValueOnly(pos, blindLeaf->getValue(pos).blind());
408  }
409  }
410  });
411  // Copy tiles serially.
412  const auto blindAcc = mBlindGrid->getConstAccessor();
413  auto iter = attTree->beginValueOn();
414  iter.setMaxDepth(AttTreeT::ValueOnIter::LEAF_DEPTH - 1);
415  for ( ; iter; ++iter) {
416  iter.modifyValue([&](AttType& v) { v = blindAcc.getValue(iter.getCoord()).blind(); });
417  }
418  } else {
419  // Here we exploit the fact that by design level sets have no active tiles.
420  // Only leaf voxels can be active.
421  sdfTree.reset(new SdfTreeT(blindTree, blindTree.background().visible(), TopologyCopy()));
422  for (BlindLeafIterT n = blindTree.cbeginLeaf(); n; ++n) {
423  const BlindLeafT& leaf = *n;
424  const openvdb::Coord xyz = leaf.origin();
425  // Get leafnodes that were allocated during topology construction!
426  SdfLeafT* sdfLeaf = sdfTree->probeLeaf(xyz);
427  AttLeafT* attLeaf = attTree->probeLeaf(xyz);
428  // Use linear offset (vs coordinate) access for better performance!
429  typename BlindLeafT::ValueOnCIter m=leaf.cbeginValueOn();
430  if (!m) {//no active values in leaf node so copy everything
431  for (openvdb::Index k = 0; k!=BlindLeafT::SIZE; ++k) {
432  const BlindType& v = leaf.getValue(k);
433  sdfLeaf->setValueOnly(k, v.visible());
434  attLeaf->setValueOnly(k, v.blind());
435  }
436  } else {//only copy active values (using flood fill for the inactive values)
437  for(; m; ++m) {
438  const openvdb::Index k = m.pos();
439  const BlindType& v = *m;
440  sdfLeaf->setValueOnly(k, v.visible());
441  attLeaf->setValueOnly(k, v.blind());
442  }
443  }
444  }
445  tools::signedFloodFill(*sdfTree);//required since we only transferred active voxels!
446  }
447 
448  if (mSdfGrid->empty()) {
449  mSdfGrid->setTree(sdfTree);
450  } else {
451  if (OutputIsMask) {
452  mSdfGrid->tree().topologyUnion(*sdfTree);
453  tools::prune(mSdfGrid->tree());
454  } else {
455  tools::csgUnion(mSdfGrid->tree(), *sdfTree, /*prune=*/true);
456  }
457  }
458 
460 }
461 
462 
463 ///////////////////////////////////////////////////////////
464 
465 
466 template<typename SdfGridT, typename AttributeT, typename InterrupterT>
467 template<typename ParticleListT, typename GridT>
468 struct ParticlesToLevelSet<SdfGridT, AttributeT, InterrupterT>::Raster
469 {
470  using DisableT = typename std::is_void<AttributeT>::type;
471  using ParticlesToLevelSetT = ParticlesToLevelSet<SdfGridT, AttributeT, InterrupterT>;
472  using SdfT = typename ParticlesToLevelSetT::SdfType; // type of signed distance values
473  using AttT = typename ParticlesToLevelSetT::AttType; // type of particle attribute
474  using ValueT = typename GridT::ValueType;
475  using AccessorT = typename GridT::Accessor;
476  using TreeT = typename GridT::TreeType;
477  using LeafNodeT = typename TreeT::LeafNodeType;
478  using PointPartitionerT = PointPartitioner<Index32, LeafNodeT::LOG2DIM>;
479 
480  static const bool
482  DoAttrXfer = !DisableT::value;
483 
484  /// @brief Main constructor
485  Raster(ParticlesToLevelSetT& parent, GridT* grid, const ParticleListT& particles)
486  : mParent(parent)
487  , mParticles(particles)
488  , mGrid(grid)
489  , mMap(*(mGrid->transform().baseMap()))
490  , mMinCount(0)
491  , mMaxCount(0)
492  , mIsCopy(false)
493  {
494  mPointPartitioner = new PointPartitionerT;
495  mPointPartitioner->construct(particles, mGrid->transform());
496  }
497 
498  /// @brief Copy constructor called by tbb threads
499  Raster(Raster& other, tbb::split)
500  : mParent(other.mParent)
501  , mParticles(other.mParticles)
502  , mGrid(new GridT(*other.mGrid, openvdb::ShallowCopy()))
503  , mMap(other.mMap)
504  , mMinCount(0)
505  , mMaxCount(0)
506  , mTask(other.mTask)
507  , mIsCopy(true)
508  , mPointPartitioner(other.mPointPartitioner)
509  {
510  mGrid->newTree();
511  }
512 
513  virtual ~Raster()
514  {
515  // Copy-constructed Rasters own temporary grids that have to be deleted,
516  // while the original has ownership of the bucket array.
517  if (mIsCopy) {
518  delete mGrid;
519  } else {
520  delete mPointPartitioner;
521  }
522  }
523 
524  void rasterizeSpheres()
525  {
526  mMinCount = mMaxCount = 0;
527  if (mParent.mInterrupter) {
528  mParent.mInterrupter->start("Rasterizing particles to level set using spheres");
529  }
530  mTask = std::bind(&Raster::rasterSpheres, std::placeholders::_1, std::placeholders::_2);
531  this->cook();
532  if (mParent.mInterrupter) mParent.mInterrupter->end();
533  }
534 
535  void rasterizeSpheres(Real radius)
536  {
537  mMinCount = radius < mParent.mRmin ? mParticles.size() : 0;
538  mMaxCount = radius > mParent.mRmax ? mParticles.size() : 0;
539  if (mMinCount>0 || mMaxCount>0) {//skipping all particles!
540  mParent.mMinCount = mMinCount;
541  mParent.mMaxCount = mMaxCount;
542  } else {
543  if (mParent.mInterrupter) {
544  mParent.mInterrupter->start(
545  "Rasterizing particles to level set using const spheres");
546  }
547  mTask = std::bind(&Raster::rasterFixedSpheres,
548  std::placeholders::_1, std::placeholders::_2, radius);
549  this->cook();
550  if (mParent.mInterrupter) mParent.mInterrupter->end();
551  }
552  }
553 
554  void rasterizeTrails(Real delta=1.0)
555  {
556  mMinCount = mMaxCount = 0;
557  if (mParent.mInterrupter) {
558  mParent.mInterrupter->start("Rasterizing particles to level set using trails");
559  }
560  mTask = std::bind(&Raster::rasterTrails,
561  std::placeholders::_1, std::placeholders::_2, delta);
562  this->cook();
563  if (mParent.mInterrupter) mParent.mInterrupter->end();
564  }
565 
566  /// @brief Kick off the optionally multithreaded computation.
567  void operator()(const tbb::blocked_range<size_t>& r)
568  {
569  assert(mTask);
570  mTask(this, r);
571  mParent.mMinCount = mMinCount;
572  mParent.mMaxCount = mMaxCount;
573  }
574 
575  /// @brief Required by tbb::parallel_reduce
576  void join(Raster& other)
577  {
579  if (OutputIsMask) {
580  if (DoAttrXfer) {
581  tools::compMax(*mGrid, *other.mGrid);
582  } else {
583  mGrid->topologyUnion(*other.mGrid);
584  }
585  } else {
586  tools::csgUnion(*mGrid, *other.mGrid, /*prune=*/true);
587  }
589  mMinCount += other.mMinCount;
590  mMaxCount += other.mMaxCount;
591  }
592 
593 private:
594  /// Disallow assignment since some of the members are references
595  Raster& operator=(const Raster&) { return *this; }
596 
597  /// @return true if the particle is too small or too large
598  bool ignoreParticle(Real R)
599  {
600  if (R < mParent.mRmin) {// below the cutoff radius
601  ++mMinCount;
602  return true;
603  }
604  if (R > mParent.mRmax) {// above the cutoff radius
605  ++mMaxCount;
606  return true;
607  }
608  return false;
609  }
610 
611  /// @brief Threaded rasterization of particles as spheres with variable radius
612  /// @param r range of indices into the list of particles
613  void rasterSpheres(const tbb::blocked_range<size_t>& r)
614  {
615  AccessorT acc = mGrid->getAccessor(); // local accessor
616  bool run = true;
617  const Real invDx = 1 / mParent.mDx;
618  AttT att;
619  Vec3R pos;
620  Real rad;
621 
622  // Loop over buckets
623  for (size_t n = r.begin(), N = r.end(); n < N; ++n) {
624  // Loop over particles in bucket n.
625  typename PointPartitionerT::IndexIterator iter = mPointPartitioner->indices(n);
626  for ( ; run && iter; ++iter) {
627  const Index32& id = *iter;
628  mParticles.getPosRad(id, pos, rad);
629  const Real R = invDx * rad;// in voxel units
630  if (this->ignoreParticle(R)) continue;
631  const Vec3R P = mMap.applyInverseMap(pos);
632  this->getAtt<DisableT>(id, att);
633  run = this->makeSphere(P, R, att, acc);
634  }//end loop over particles
635  }//end loop over buckets
636  }
637 
638  /// @brief Threaded rasterization of particles as spheres with a fixed radius
639  /// @param r range of indices into the list of particles
640  /// @param R radius of fixed-size spheres
641  void rasterFixedSpheres(const tbb::blocked_range<size_t>& r, Real R)
642  {
643  AccessorT acc = mGrid->getAccessor(); // local accessor
644  AttT att;
645  Vec3R pos;
646 
647  // Loop over buckets
648  for (size_t n = r.begin(), N = r.end(); n < N; ++n) {
649  // Loop over particles in bucket n.
650  for (auto iter = mPointPartitioner->indices(n); iter; ++iter) {
651  const Index32& id = *iter;
652  this->getAtt<DisableT>(id, att);
653  mParticles.getPos(id, pos);
654  const Vec3R P = mMap.applyInverseMap(pos);
655  this->makeSphere(P, R, att, acc);
656  }
657  }
658  }
659 
660  /// @brief Threaded rasterization of particles as spheres with velocity trails
661  /// @param r range of indices into the list of particles
662  /// @param delta inter-sphere spacing
663  void rasterTrails(const tbb::blocked_range<size_t>& r, Real delta)
664  {
665  AccessorT acc = mGrid->getAccessor(); // local accessor
666  bool run = true;
667  AttT att;
668  Vec3R pos, vel;
669  Real rad;
670  const Vec3R origin = mMap.applyInverseMap(Vec3R(0,0,0));
671  const Real Rmin = mParent.mRmin, invDx = 1 / mParent.mDx;
672 
673  // Loop over buckets
674  for (size_t n = r.begin(), N = r.end(); n < N; ++n) {
675  // Loop over particles in bucket n.
676  typename PointPartitionerT::IndexIterator iter = mPointPartitioner->indices(n);
677  for ( ; run && iter; ++iter) {
678  const Index32& id = *iter;
679  mParticles.getPosRadVel(id, pos, rad, vel);
680  const Real R0 = invDx * rad;
681  if (this->ignoreParticle(R0)) continue;
682  this->getAtt<DisableT>(id, att);
683  const Vec3R P0 = mMap.applyInverseMap(pos);
684  const Vec3R V = mMap.applyInverseMap(vel) - origin; // exclude translation
685  const Real speed = V.length(), invSpeed = 1.0 / speed;
686  const Vec3R Nrml = -V * invSpeed; // inverse normalized direction
687  Vec3R P = P0; // local position of instance
688  Real R = R0, d = 0; // local radius and length of trail
689  while (run && d <= speed) {
690  run = this->makeSphere(P, R, att, acc);
691  P += 0.5 * delta * R * Nrml; // adaptive offset along inverse velocity direction
692  d = (P - P0).length(); // current length of trail
693  R = R0 - (R0 - Rmin) * d * invSpeed; // R = R0 -> mRmin(e.g. 1.5)
694  }//end loop over sphere instances
695  }//end loop over particles
696  }//end loop over buckets
697  }
698 
699  void cook()
700  {
701  // parallelize over the point buckets
702  const Index32 bucketCount = Index32(mPointPartitioner->size());
703 
704  if (mParent.mGrainSize>0) {
705  tbb::parallel_reduce(
706  tbb::blocked_range<size_t>(0, bucketCount, mParent.mGrainSize), *this);
707  } else {
708  (*this)(tbb::blocked_range<size_t>(0, bucketCount));
709  }
710  }
711 
712  /// @brief Rasterize sphere at position P and radius R into
713  /// a narrow-band level set with half-width, mHalfWidth.
714  /// @return @c false if rasterization was interrupted
715  ///
716  /// @param P coordinates of the particle position in voxel units
717  /// @param R radius of particle in voxel units
718  /// @param att an optional user-defined attribute value to be associated with voxels
719  /// @param acc grid accessor with a private copy of the grid
720  ///
721  /// @note For best performance all computations are performed in voxel space,
722  /// with the important exception of the final level set value that is converted
723  /// to world units (the grid stores the closest Euclidean signed distances
724  /// measured in world units). Also note we use the convention of positive distances
725  /// outside the surface and negative distances inside the surface.
726  template <bool IsMaskT = OutputIsMask>
728  makeSphere(const Vec3R& P, Real R, const AttT& att, AccessorT& acc)
729  {
730  const Real
731  dx = mParent.mDx,
732  w = mParent.mHalfWidth,
733  max = R + w, // maximum distance in voxel units
734  max2 = math::Pow2(max), // square of maximum distance in voxel units
735  min2 = math::Pow2(math::Max(Real(0), R - w)); // square of minimum distance
736  // Bounding box of the sphere
737  const Coord
738  lo(math::Floor(P[0]-max),math::Floor(P[1]-max),math::Floor(P[2]-max)),
739  hi(math::Ceil( P[0]+max),math::Ceil( P[1]+max),math::Ceil( P[2]+max));
740  const ValueT inside = -mGrid->background();
741 
742  ValueT v;
743  size_t count = 0;
744  for (Coord c = lo; c.x() <= hi.x(); ++c.x()) {
745  //only check interrupter every 32'th scan in x
746  if (!(count++ & ((1<<5)-1)) && util::wasInterrupted(mParent.mInterrupter)) {
747  thread::cancelGroupExecution();
748  return false;
749  }
750  const Real x2 = math::Pow2(c.x() - P[0]);
751  for (c.y() = lo.y(); c.y() <= hi.y(); ++c.y()) {
752  const Real x2y2 = x2 + math::Pow2(c.y() - P[1]);
753  for (c.z() = lo.z(); c.z() <= hi.z(); ++c.z()) {
754  const Real x2y2z2 = x2y2 + math::Pow2(c.z()-P[2]); // squared dist from c to P
755 #if defined __INTEL_COMPILER
756  _Pragma("warning (push)")
757  _Pragma("warning (disable:186)") // "pointless comparison of unsigned integer with zero"
758 #endif
759  if (x2y2z2 >= max2 || (!acc.probeValue(c, v) && (v < ValueT(0))))
760  continue;//outside narrow band of the particle or inside existing level set
761 #if defined __INTEL_COMPILER
762  _Pragma("warning (pop)")
763 #endif
764  if (x2y2z2 <= min2) {//inside narrow band of the particle.
765  acc.setValueOff(c, inside);
766  continue;
767  }
768  // convert signed distance from voxel units to world units
769  //const ValueT d=dx*(math::Sqrt(x2y2z2) - R);
770  const ValueT d = Merge(static_cast<SdfT>(dx*(math::Sqrt(x2y2z2)-R)), att);
771  if (d < v) acc.setValue(c, d);//CSG union
772  }//end loop over z
773  }//end loop over y
774  }//end loop over x
775  return true;
776  }
777 
778  /// @brief Rasterize a sphere of radius @a r at position @a p into a boolean mask grid.
779  /// @return @c false if rasterization was interrupted
780  template <bool IsMaskT = OutputIsMask>
782  makeSphere(const Vec3R& p, Real r, const AttT& att, AccessorT& acc)
783  {
784  const Real
785  rSquared = r * r, // sphere radius squared, in voxel units
786  inW = r / math::Sqrt(6.0); // half the width in voxel units of an inscribed cube
787  const Coord
788  // Bounding box of the sphere
789  outLo(math::Floor(p[0] - r), math::Floor(p[1] - r), math::Floor(p[2] - r)),
790  outHi(math::Ceil(p[0] + r), math::Ceil(p[1] + r), math::Ceil(p[2] + r)),
791  // Bounds of the inscribed cube
792  inLo(math::Ceil(p[0] - inW), math::Ceil(p[1] - inW), math::Ceil(p[2] - inW)),
793  inHi(math::Floor(p[0] + inW), math::Floor(p[1] + inW), math::Floor(p[2] + inW));
794  // Bounding boxes of regions comprising out - in
795  /// @todo These could be divided further into sparsely- and densely-filled subregions.
796  const std::vector<CoordBBox> padding{
797  CoordBBox(outLo.x(), outLo.y(), outLo.z(), inLo.x()-1, outHi.y(), outHi.z()),
798  CoordBBox(inHi.x()+1, outLo.y(), outLo.z(), outHi.x(), outHi.y(), outHi.z()),
799  CoordBBox(outLo.x(), outLo.y(), outLo.z(), outHi.x(), inLo.y()-1, outHi.z()),
800  CoordBBox(outLo.x(), inHi.y()+1, outLo.z(), outHi.x(), outHi.y(), outHi.z()),
801  CoordBBox(outLo.x(), outLo.y(), outLo.z(), outHi.x(), outHi.y(), inLo.z()-1),
802  CoordBBox(outLo.x(), outLo.y(), inHi.z()+1, outHi.x(), outHi.y(), outHi.z()),
803  };
804  const ValueT onValue = Merge(SdfT(1), att);
805 
806  // Sparsely fill the inscribed cube.
807  /// @todo Use sparse fill only if 2r > leaf width?
808  acc.tree().sparseFill(CoordBBox(inLo, inHi), onValue);
809 
810  // Densely fill the remaining regions.
811  for (const auto& bbox: padding) {
812  if (util::wasInterrupted(mParent.mInterrupter)) {
813  thread::cancelGroupExecution();
814  return false;
815  }
816  const Coord &bmin = bbox.min(), &bmax = bbox.max();
817  Coord c;
818  Real cx, cy, cz;
819  for (c = bmin, cx = c.x(); c.x() <= bmax.x(); ++c.x(), cx += 1) {
820  const Real x2 = math::Pow2(cx - p[0]);
821  for (c.y() = bmin.y(), cy = c.y(); c.y() <= bmax.y(); ++c.y(), cy += 1) {
822  const Real x2y2 = x2 + math::Pow2(cy - p[1]);
823  for (c.z() = bmin.z(), cz = c.z(); c.z() <= bmax.z(); ++c.z(), cz += 1) {
824  const Real x2y2z2 = x2y2 + math::Pow2(cz - p[2]);
825  if (x2y2z2 < rSquared) {
826  acc.setValue(c, onValue);
827  }
828  }
829  }
830  }
831  }
832  return true;
833  }
834 
835  using FuncType = typename std::function<void (Raster*, const tbb::blocked_range<size_t>&)>;
836 
837  template<typename DisableType>
839  getAtt(size_t, AttT&) const {}
840 
841  template<typename DisableType>
843  getAtt(size_t n, AttT& a) const { mParticles.getAtt(n, a); }
844 
845  template<typename T>
847  Merge(T s, const AttT&) const { return s; }
848 
849  template<typename T>
851  Merge(T s, const AttT& a) const { return ValueT(s,a); }
852 
853  ParticlesToLevelSetT& mParent;
854  const ParticleListT& mParticles;//list of particles
855  GridT* mGrid;
856  const math::MapBase& mMap;
857  size_t mMinCount, mMaxCount;//counters for ignored particles!
858  FuncType mTask;
859  const bool mIsCopy;
860  PointPartitionerT* mPointPartitioner;
861 }; // struct ParticlesToLevelSet::Raster
862 
863 
864 ///////////////////// YOU CAN SAFELY IGNORE THIS SECTION /////////////////////
865 
866 /// @cond OPENVDB_DOCS_INTERNAL
867 
868 namespace p2ls_internal {
869 
870 // This is a simple type that combines a distance value and a particle
871 // attribute. It's required for attribute transfer which is defined in the
872 // Raster class above.
873 /// @private
874 template<typename VisibleT, typename BlindT>
875 class BlindData
876 {
877 public:
878  using type = VisibleT;
879  using VisibleType = VisibleT;
880  using BlindType = BlindT;
881 
882  BlindData() {}
883  explicit BlindData(VisibleT v) : mVisible(v), mBlind(zeroVal<BlindType>()) {}
884  BlindData(VisibleT v, BlindT b) : mVisible(v), mBlind(b) {}
885  BlindData(const BlindData&) = default;
886  BlindData& operator=(const BlindData&) = default;
887  const VisibleT& visible() const { return mVisible; }
888  const BlindT& blind() const { return mBlind; }
890  bool operator==(const BlindData& rhs) const { return mVisible == rhs.mVisible; }
892  bool operator< (const BlindData& rhs) const { return mVisible < rhs.mVisible; }
893  bool operator> (const BlindData& rhs) const { return mVisible > rhs.mVisible; }
894  BlindData operator+(const BlindData& rhs) const { return BlindData(mVisible + rhs.mVisible); }
895  BlindData operator-(const BlindData& rhs) const { return BlindData(mVisible - rhs.mVisible); }
896  BlindData operator-() const { return BlindData(-mVisible, mBlind); }
897 
898 protected:
899  VisibleT mVisible;
900  BlindT mBlind;
901 };
902 
903 /// @private
904 // Required by several of the tree nodes
905 template<typename VisibleT, typename BlindT>
906 inline std::ostream& operator<<(std::ostream& ostr, const BlindData<VisibleT, BlindT>& rhs)
907 {
908  ostr << rhs.visible();
909  return ostr;
910 }
911 
912 /// @private
913 // Required by math::Abs
914 template<typename VisibleT, typename BlindT>
915 inline BlindData<VisibleT, BlindT> Abs(const BlindData<VisibleT, BlindT>& x)
916 {
917  return BlindData<VisibleT, BlindT>(math::Abs(x.visible()), x.blind());
918 }
919 
920 /// @private
921 // Required to support the (zeroVal<BlindData>() + val) idiom.
922 template<typename VisibleT, typename BlindT, typename T>
923 inline BlindData<VisibleT, BlindT>
924 operator+(const BlindData<VisibleT, BlindT>& x, const T& rhs)
925 {
926  return BlindData<VisibleT, BlindT>(x.visible() + static_cast<VisibleT>(rhs), x.blind());
927 }
928 
929 } // namespace p2ls_internal
930 
931 /// @endcond
932 
933 //////////////////////////////////////////////////////////////////////////////
934 
935 
936 // The following are convenience functions for common use cases.
937 
938 template<typename GridT, typename ParticleListT, typename InterrupterT>
939 inline void
940 particlesToSdf(const ParticleListT& plist, GridT& grid, InterrupterT* interrupt)
941 {
943  "particlesToSdf requires an SDF grid with floating-point values");
944 
945  if (grid.getGridClass() != GRID_LEVEL_SET) {
946  OPENVDB_LOG_WARN("particlesToSdf requires a level set grid;"
947  " try Grid::setGridClass(openvdb::GRID_LEVEL_SET)");
948  }
949 
950  ParticlesToLevelSet<GridT> p2ls(grid, interrupt);
951  p2ls.rasterizeSpheres(plist);
952  tools::pruneLevelSet(grid.tree());
953 }
954 
955 template<typename GridT, typename ParticleListT, typename InterrupterT>
956 inline void
957 particlesToSdf(const ParticleListT& plist, GridT& grid, Real radius, InterrupterT* interrupt)
958 {
960  "particlesToSdf requires an SDF grid with floating-point values");
961 
962  if (grid.getGridClass() != GRID_LEVEL_SET) {
963  OPENVDB_LOG_WARN("particlesToSdf requires a level set grid;"
964  " try Grid::setGridClass(openvdb::GRID_LEVEL_SET)");
965  }
966 
967  ParticlesToLevelSet<GridT> p2ls(grid, interrupt);
968  p2ls.rasterizeSpheres(plist, radius);
969  tools::pruneLevelSet(grid.tree());
970 }
971 
972 template<typename GridT, typename ParticleListT, typename InterrupterT>
973 inline void
974 particleTrailsToSdf(const ParticleListT& plist, GridT& grid, Real delta, InterrupterT* interrupt)
975 {
977  "particleTrailsToSdf requires an SDF grid with floating-point values");
978 
979  if (grid.getGridClass() != GRID_LEVEL_SET) {
980  OPENVDB_LOG_WARN("particlesToSdf requires a level set grid;"
981  " try Grid::setGridClass(openvdb::GRID_LEVEL_SET)");
982  }
983 
984  ParticlesToLevelSet<GridT> p2ls(grid, interrupt);
985  p2ls.rasterizeTrails(plist, delta);
986  tools::pruneLevelSet(grid.tree());
987 }
988 
989 template<typename GridT, typename ParticleListT, typename InterrupterT>
990 inline void
991 particlesToMask(const ParticleListT& plist, GridT& grid, InterrupterT* interrupt)
992 {
994  "particlesToMask requires a boolean-valued grid");
995  ParticlesToLevelSet<GridT> p2ls(grid, interrupt);
996  p2ls.rasterizeSpheres(plist);
997  tools::prune(grid.tree());
998 }
999 
1000 template<typename GridT, typename ParticleListT, typename InterrupterT>
1001 inline void
1002 particlesToMask(const ParticleListT& plist, GridT& grid, Real radius, InterrupterT* interrupt)
1003 {
1005  "particlesToMask requires a boolean-valued grid");
1006  ParticlesToLevelSet<GridT> p2ls(grid, interrupt);
1007  p2ls.rasterizeSpheres(plist, radius);
1008  tools::prune(grid.tree());
1009 }
1010 
1011 template<typename GridT, typename ParticleListT, typename InterrupterT>
1012 inline void
1013 particleTrailsToMask(const ParticleListT& plist, GridT& grid, Real delta, InterrupterT* interrupt)
1014 {
1016  "particleTrailsToMask requires a boolean-valued grid");
1017  ParticlesToLevelSet<GridT> p2ls(grid, interrupt);
1018  p2ls.rasterizeTrails(plist, delta);
1019  tools::prune(grid.tree());
1020 }
1021 
1022 } // namespace tools
1023 } // namespace OPENVDB_VERSION_NAME
1024 } // namespace openvdb
1025 
1026 #endif // OPENVDB_TOOLS_PARTICLES_TO_LEVELSET_HAS_BEEN_INCLUDED
size_t getMaxCount() const
Return the number of particles that were ignored because they were larger than the maximum radius...
Real getHalfWidth() const
Return the half-width of the narrow band in voxel units.
int Ceil(float x)
Return the ceiling of x.
Definition: Math.h:856
MeshToVoxelEdgeData::EdgeData Abs(const MeshToVoxelEdgeData::EdgeData &x)
Functions to efficiently perform various compositing operations on grids.
math::Vec3< Real > Vec3R
Definition: Types.h:72
Tag dispatch class that distinguishes shallow copy constructors from deep copy constructors.
Definition: Types.h:680
Type Pow2(Type x)
Return x2.
Definition: Math.h:548
Real getRmax() const
Return the largest radius allowed in voxel units.
void
Definition: png.h:1083
AttGridType::Ptr attributeGrid()
Return a pointer to the grid containing the optional user-defined attribute.
Real getVoxelSize() const
Return the size of a voxel in world units.
const GLdouble * v
Definition: glcorearb.h:837
GLsizei const GLfloat * value
Definition: glcorearb.h:824
void finalize(bool prune=false)
This method syncs up the level set and attribute grids and therefore needs to be called before any of...
void particlesToSdf(const ParticleListT &, GridT &, InterrupterT *=nullptr)
Populate a scalar, floating-point grid with CSG-unioned level set spheres described by the given part...
IMATH_HOSTDEVICE constexpr Plane3< T > operator-(const Plane3< T > &plane) IMATH_NOEXCEPT
Reflect the pla.
Definition: ImathPlane.h:253
#define OPENVDB_LOG_WARN(mesg)
Definition: logging.h:277
GLboolean GLboolean GLboolean GLboolean a
Definition: glcorearb.h:1222
GLdouble s
Definition: glad.h:3009
void csgUnion(GridOrTreeT &a, GridOrTreeT &b, bool prune=true, bool pruneCancelledTiles=false)
Given two level set grids, replace the A grid with the union of A and B.
Definition: Composite.h:886
GLuint GLsizei GLsizei * length
Definition: glcorearb.h:795
#define OPENVDB_USE_VERSION_NAMESPACE
Definition: version.h:239
void particlesToMask(const ParticleListT &, GridT &, InterrupterT *=nullptr)
Activate a boolean grid wherever it intersects the spheres described by the given particle positions ...
GA_API const UT_StringHolder P
void setGrainSize(int grainSize)
Set the grain size used for threading.
GLdouble GLdouble x2
Definition: glad.h:2349
ParticlesToLevelSet(SdfGridT &grid, InterrupterT *interrupt=nullptr)
Constructor using an existing boolean or narrow-band level set grid.
Spatially partitions points using a parallel radix-based sorting algorithm.
SYS_FORCE_INLINE const_iterator end() const
bool operator==(const BaseDimensions< T > &a, const BaseDimensions< Y > &b)
Definition: Dimensions.h:137
constexpr T zeroVal()
Return the value of type T that corresponds to zero.
Definition: Math.h:70
int getGrainSize() const
Return the grain size used for threading.
OIIO_FORCEINLINE vbool4 operator>(const vint4 &a, const vint4 &b)
Definition: simd.h:4561
GLdouble n
Definition: glcorearb.h:2008
Coord Abs(const Coord &xyz)
Definition: Coord.h:517
void rasterizeSpheres(const ParticleListT &pa)
Rasterize each particle as a sphere with the particle's position and radius.
typename std::conditional< DisableT::value, size_t, AttributeT >::type AttType
General-purpose arithmetic and comparison routines, most of which accept arbitrary value types (or at...
void setRmax(Real Rmax)
Set the largest radius allowed in voxel units.
float Sqrt(float x)
Return the square root of a floating-point value.
Definition: Math.h:761
Defined various multi-threaded utility functions for trees.
This class manages a linear array of pointers to a given tree's leaf nodes, as well as optional auxil...
Definition: LeafManager.h:84
bool operator<(const GU_TetrahedronFacet &a, const GU_TetrahedronFacet &b)
GLuint id
Definition: glcorearb.h:655
#define OPENVDB_NO_FP_EQUALITY_WARNING_BEGIN
Definition: Math.h:48
void foreach(const LeafOp &op, bool threaded=true, size_t grainSize=1)
Threaded method that applies a user-supplied functor to each leaf node in the LeafManager.
Definition: LeafManager.h:483
typename SdfGridT::template ValueConverter< AttType >::Type AttGridType
GLboolean GLboolean GLboolean b
Definition: glcorearb.h:1222
GA_API const UT_StringHolder transform
GLint GLenum GLint x
Definition: glcorearb.h:409
IMATH_HOSTDEVICE constexpr Quat< T > operator+(const Quat< T > &q1, const Quat< T > &q2) IMATH_NOEXCEPT
Quaterion addition.
Definition: ImathQuat.h:905
size_t getMinCount() const
Return the number of particles that were ignored because they were smaller than the minimum radius...
#define OPENVDB_NO_UNREACHABLE_CODE_WARNING_BEGIN
Definition: Platform.h:146
Propagate the signs of distance values from the active voxels in the narrow band to the inactive valu...
bool ignoredParticles() const
Return true if any particles were ignored due to their size.
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 ...
Definition: Composite.h:748
void particleTrailsToSdf(const ParticleListT &, GridT &, Real delta=1, InterrupterT *=nullptr)
Populate a scalar, floating-point grid with CSG-unioned trails of level set spheres with decreasing r...
void setRmin(Real Rmin)
Set the smallest radius allowed in voxel units.
void rasterizeTrails(const ParticleListT &pa, Real delta=1.0)
Rasterize each particle as a trail comprising the CSG union of spheres of decreasing radius...
void signedFloodFill(TreeOrLeafManagerT &tree, bool threaded=true, size_t grainSize=1, Index minLevel=0)
Set the values of all inactive voxels and tiles of a narrow-band level set from the signs of the acti...
ImageBuf OIIO_API max(Image_or_Const A, Image_or_Const B, ROI roi={}, int nthreads=0)
GA_API const UT_StringHolder N
void particleTrailsToMask(const ParticleListT &, GridT &, Real delta=1, InterrupterT *=nullptr)
Activate a boolean grid wherever it intersects trails of spheres with decreasing radius, where the starting position and radius and the direction of each trail is given by particle attributes.
#define SIZE
Definition: simple.C:41
#define OPENVDB_NO_UNREACHABLE_CODE_WARNING_END
Definition: Platform.h:147
int Floor(float x)
Return the floor of x.
Definition: Math.h:848
GLubyte GLubyte GLubyte GLubyte w
Definition: glcorearb.h:857
Tag dispatch class that distinguishes topology copy constructors from deep copy constructors.
Definition: Types.h:683
A LeafManager manages a linear array of pointers to a given tree's leaf nodes, as well as optional au...
GLboolean r
Definition: glcorearb.h:1222
void prune(TreeT &tree, typename TreeT::ValueType tolerance=zeroVal< typename TreeT::ValueType >(), bool threaded=true, size_t grainSize=1)
Reduce the memory footprint of a tree by replacing with tiles any nodes whose values are all the same...
Definition: Prune.h:335
void OIIO_UTIL_API split(string_view str, std::vector< string_view > &result, string_view sep=string_view(), int maxsplit=-1)
void pruneLevelSet(TreeT &tree, bool threaded=true, size_t grainSize=1)
Reduce the memory footprint of a tree by replacing nodes whose values are all inactive with inactive ...
Definition: Prune.h:390
typename std::is_void< AttributeT >::type DisableT
bool wasInterrupted(T *i, int percent=-1)
type
Definition: core.h:1059
#define OPENVDB_VERSION_NAME
The version namespace name for this library version.
Definition: version.h:119
Real getRmin() const
Return the smallest radius allowed in voxel units.
const Type & Max(const Type &a, const Type &b)
Return the maximum of two values.
Definition: Math.h:595
GLint GLsizei count
Definition: glcorearb.h:405
#define OPENVDB_THROW(exception, message)
Definition: Exceptions.h:74
auto join(It begin, Sentinel end, string_view sep) -> join_view< It, Sentinel >
Definition: format.h:2559
#define OPENVDB_NO_FP_EQUALITY_WARNING_END
Definition: Math.h:49