VTK  9.6.20260213
vtkVoronoiCore3D.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
100
101#ifndef vtkVoronoiCore3D_h
102#define vtkVoronoiCore3D_h
103
104#include "vtkLocatorInterface.h" // Interface to spatial locators
105#include "vtkPoints.h" // for input points
106#include "vtkSMPTools.h" // SMP parallel processing
107#include "vtkShellBinIterator.h" // Iteration over the bins of the locator
108#include "vtkStaticPointLocator.h" // for locating nearby points
109#include "vtkVoronoiCore.h" // Core includes
110#include "vtkVoronoiHull.h" // Generate Voronoi convex hulls
111
112VTK_ABI_NAMESPACE_BEGIN
113
128
141{
142 // Optional regions ids for point classification.
143 const int* Regions;
144
145 // Constructors
147 : Regions(nullptr)
148 {
149 }
150 vtkVoronoiClassifier3D(const int* regions)
151 : Regions(regions)
152 {
153 }
154
155 // Method required by vtkVoronoiCore3D.
157 {
158 if (c)
159 {
160 this->Regions = c->Regions;
161 }
162 }
163
164 // Method required by vtkVoronoiCore3D.
166 vtkVoronoiSpokesType& spokes, int& numSpokes, int& maxPoints, int& maxFaces);
167
168 // Method required by vtkVoronoiCore3D. By default, any region id >= 0 is
169 // considered a valid inside region (<0 region values are reserved for
170 // algorithm use). If no region ids have been specified, than the point
171 // ptId is inside an interior region.
173 {
174 if (!this->Regions)
175 {
176 return (ptId >= 0);
177 }
178 else
179 {
180 return (ptId >= 0 && this->Regions[ptId] >= 0);
181 }
182 }
183
184 // Method required by vtkVoronoiCore3D. Determine if the two points
185 // ptId and neiId (which form a spoke) are in the same region. It is
186 // assumed that both (ptId,neiId >= 0), i.e., inside.
188 {
189 return (!this->Regions || this->Regions[ptId] == this->Regions[neiId]);
190 }
191}; // vtkVoronoiClassifier3D
192
197template <class TCompositor, class TClassifier>
199{
200 vtkIdType ThreadId; // assign a thread id [0,NumThreadsUsed)
201 int MaxPoints; // the maximum number of points in any hull
202 int MaxFaces; // the maximum number of faces in any hull
203 int NumPrunes; // total number of pruning operations
204 vtkBatchIdsType LocalBatches; // list of batches processed by this thread
205 vtkVoronoiSpokesType LocalSpokes; // connecting edges/spokes for each hull
206 vtkShellBinIterator SIter; // iterator over static point locator bins
207 vtkVoronoiHull Hull; // computational 3D Voronoi hull algorithm
208 typename TCompositor::LocalData Compositor; // gather data from compositing operations
209 TClassifier Classifier; // Used to classify spokes (based on regions)
210
212 : ThreadId(-1)
213 , MaxPoints(0)
214 , MaxFaces(0)
215 , NumPrunes(0)
216 {
217 this->LocalBatches.reserve(2048);
218 this->LocalSpokes.reserve(2048);
219 }
220}; // vtkVoronoi3DLocalData
221
226template <class TCompositor, class TClassifier>
227using ThreadMapType = std::vector<vtkVoronoi3DLocalData<TCompositor, TClassifier>*>;
228
241template <class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
243{
244public:
261 static std::unique_ptr<vtkVoronoiCore3D<TCompositor, TClassifier>> Execute(vtkAlgorithm* filter,
262 unsigned int batchSize, vtkStaticPointLocator* loc, vtkPoints* inPts, double padding,
263 vtkIdType maxClips, bool validate, double pruneTol, TCompositor* comp, TClassifier* cl);
264
266
271 int GetNumberOfThreads() { return this->NumberOfThreads; }
273 {
274 return this->ThreadMap[threadNum];
275 }
276
277
279
285 int GetMaximumNumberOfPoints() { return this->MaximumNumberOfPoints; }
286 int GetMaximumNumberOfFaces() { return this->MaximumNumberOfFaces; }
287 int GetNumberOfPrunes() { return this->NumberOfPrunes; }
289
295
301 vtkStaticPointLocator* loc, vtkPoints* inPts, double padding, vtkIdType maxClips, bool validate,
302 double pruneTol, TCompositor* comp, TClassifier* c);
303
305
313 void operator()(vtkIdType batchId, vtkIdType endBatchId);
314 void Reduce();
316
324 TCompositor Compositor;
325
331 TClassifier Classifier;
332
339
346
351 vtkIdType GetNumberOfPoints() { return this->NPts; }
352 const double* GetPoints() { return this->Points; }
353
359
360private:
361 vtkIdType NPts; // The number of input points
362 vtkPoints* InPoints; // Input points as data array
363 const double* Points; // Input points as pointer to x-y-z AOS doubles
364 vtkStaticPointLocator* Locator; // Used to (quickly) find nearby points
365 double Padding; // The padding distance around the bounding box
366 double Bounds[6]; // locator bounds
367 double PaddedBounds[6]; // the domain over which Voronoi is calculated
368 vtkIdType MaxClips; // Limit the number of half-space clips
369
370 // Enable pruning of spokes (equivalent to deletion of a degenerate hull face)
371 bool Validate; // Indicate whether to explicitly validate the mesh
372 vtkIdType NumberOfPrunes; // If pruning is on, keep track of the number of prunes
373 double PruneTolerance; // Specify the spoke prune tolerance
374
375 // High-level information captured during processing
376 int NumberOfThreads; // Keep track of the number of threads used durinf processing
377 int MaximumNumberOfPoints; // Maximum numper of points in a generated Voronoi hull
378 int MaximumNumberOfFaces; // Maximum number of faces (i.e., spokes) in a generated Voronoi hull
379
380 // Storage local to each thread, as well as working/scratch arrays. We
381 // don't want to allocate working arrays on every thread invocation. Thread
382 // local storage saves lots of new/delete (e.g. the locator tuples).
385
390 using vtkBinIterator = vtkShellBinIterator;
391 bool BuildHull(vtkVoronoiHull& hull, vtkBinIterator* siter, const double* pts, vtkIdType maxClips,
392 vtkDist2TupleArray& results, int& numPrunes);
393
394public:
404 {
407 void operator()(vtkIdType threadId, vtkIdType endThreadId);
408
409 // Invoke the production of wheels and spokes
411 }; // ProduceWheelsAndSpokes
412
426 {
429
430 vtkMergeTuples3DType MergeTuples; // temporary array for merging points
431 vtkMergeMapType MergeMap; // maps tile/hull point ids to merged point ids
432 vtkIdType NumMergedPts; // after merging, the number of points remaining
433
440 vtkIdType GetNumberOfMergedPoints() { return this->NumMergedPts; }
441
442 // Core Voronoi methods in support of vtkSMPTools
443 void Initialize() {}
444 void operator()(vtkIdType threadId, vtkIdType endThreadId);
445 void Reduce();
446
447 // Execute the topological point merge to produce a merge map.
448 static std::unique_ptr<TopologicalMerge> Execute(
450 }; // TopologicalMerge
451
452}; // vtkVoronoiCore3D
453
458
467{
473 void Finalize() {}
474
479 {
481 void AddData(vtkVoronoiHull& vtkNotUsed(hull), int vtkNotUsed(numSpokes),
482 const vtkVoronoiSpoke* vtkNotUsed(spokes))
483 {
484 }
485 };
486}; // vtkEmptyVoronoi3DCompositor
487
488// Almost minimal classifier - just records the hull's number of points and faces. It also
489// considers regions if any are defined -- this is optional in some cases.
491{
492 // Optional region ids
493 const int* Regions;
494
495 // Constructor
497 : Regions(nullptr)
498 {
499 }
500 vtkEmptyVoronoi3DClassifier(const int* regions)
501 : Regions(regions)
502 {
503 }
504
505 // Initialize
507 {
508 if (c)
509 {
510 this->Regions = c->Regions;
511 }
512 }
513
514 // Method required by vtkVoronoiCore3D. This vtkEmptyClassifier provides
515 // the minimum information needed.
517 vtkVoronoiSpokesType& vtkNotUsed(spokes), int& numSpokes, int& maxPoints, int& maxFaces)
518 {
519 wheels[hull.PtId] = hull.NumFaces; // numFaces == numSpokes
520 maxPoints = (hull.NumPts > maxPoints ? hull.NumPts : maxPoints);
521 maxFaces = (hull.NumFaces > maxFaces ? hull.NumFaces : maxFaces);
522 numSpokes = 0;
523 return nullptr;
524 }
525
526 // If no region ids are provided, all processed points are inside the same
527 // region.
529 {
530 if (!this->Regions)
531 {
532 return (ptId >= 0);
533 }
534 else
535 {
536 return (ptId >= 0 && this->Regions[ptId] >= 0);
537 }
538 }
540 {
541 return (!this->Regions || this->Regions[ptId] == this->Regions[neiId]);
542 }
543}; // vtkEmptyVoronoi3DClassifier
544
545VTK_ABI_NAMESPACE_END
546#include "vtkVoronoiCore3D.txx"
547
548#endif
549// VTK-HeaderTest-Exclude: vtkVoronoiCore3D.h
represent and manipulate 3D points
Definition vtkPoints.h:139
Thread local storage for VTK objects.
A fast, lightweight class for iterating over the bins of a 3D vtkStaticPointLocator.
quickly locate points in 3-space
int GetNumberOfThreads()
Access the local thread data produced by execution of the filter.
TCompositor Compositor
The compositor enables this vtkVoronoiCore3D templated class to be used in different applications.
ThreadMapType< TCompositor, TClassifier > ThreadMap
vtkVoronoiAdjacencyGraph & GetAdjacencyGraph()
Obtain the adjacency graph (wheel & spokes data structure).
TClassifier Classifier
This templated class is used to extend the API of this vtkVoronoiCore3D class, to implement the spoke...
int GetNumberOfPrunes()
Obtain information about the execution of the Voronoi algorithm.
void Initialize()
Core vtkSMPTools methods to process hulls in parallel, and produce thread local output data.
void operator()(vtkIdType batchId, vtkIdType endBatchId)
Core vtkSMPTools methods to process hulls in parallel, and produce thread local output data.
vtkVoronoiCore3D(vtkAlgorithm *filter, vtkVoronoiBatchManager &batcher, vtkStaticPointLocator *loc, vtkPoints *inPts, double padding, vtkIdType maxClips, bool validate, double pruneTol, TCompositor *comp, TClassifier *c)
Constructor.
int GetMaximumNumberOfPoints()
Obtain information about the execution of the Voronoi algorithm.
static std::unique_ptr< vtkVoronoiCore3D< TCompositor, TClassifier > > Execute(vtkAlgorithm *filter, unsigned int batchSize, vtkStaticPointLocator *loc, vtkPoints *inPts, double padding, vtkIdType maxClips, bool validate, double pruneTol, TCompositor *comp, TClassifier *cl)
A factory method to conveniently instantiate and execute the algorithm.
vtkAlgorithm * Filter
Used for controlling filter abort and accessing filter information.
vtkVoronoiAdjacencyGraph Graph
This is used to create the spokes and wheels adjacency graph used to validate the tessellation and pr...
vtkIdType GetNumberOfPoints()
Convenience methods to retrieve the number of input points, and the raw double* points array.
void Reduce()
Core vtkSMPTools methods to process hulls in parallel, and produce thread local output data.
int GetMaximumNumberOfFaces()
Obtain information about the execution of the Voronoi algorithm.
vtkVoronoi3DLocalData< TCompositor, TClassifier > * GetThreadData(int threadNum)
Access the local thread data produced by execution of the filter.
const double * GetPoints()
vtkVoronoiBatchManager Batcher
Controls processing of batches of generating points.
The polyhedron class proper.
Represent an array of vtkDist2Tuples.
bool IsInsideRegion(vtkIdType ptId)
const vtkVoronoiSpoke * AddAdjacencyInformation(vtkVoronoiHull &hull, vtkVoronoiWheelsType &wheels, vtkVoronoiSpokesType &spokes, int &numSpokes, int &maxPoints, int &maxFaces)
vtkEmptyVoronoi3DClassifier(const int *regions)
void Initialize(vtkEmptyVoronoi3DClassifier *c)
bool IsSameRegion(vtkIdType ptId, vtkIdType neiId)
Thread local data may be needed.
void AddData(vtkVoronoiHull &hull, int numSpokes, const vtkVoronoiSpoke *spokes)
void Initialize(vtkEmptyVoronoi3DCompositor *c)
These are convenience/demonstration classes for configuring the templated 3D Voronoi classes.
void Initialize(vtkIdType numPts, vtkEmptyVoronoi3DCompositor *)
Prepare to accumulate compositing information: specify the total number of generating points to be pr...
The following thread local data is used to process and keep track of information on a per-thread basi...
vtkVoronoiSpokesType LocalSpokes
TCompositor::LocalData Compositor
vtkBatchIdsType LocalBatches
vtkShellBinIterator SIter
The adjacency graph, a collection of wheels and spokes, is a topological construct that connects Voro...
Class to manage batches of points.
vtkVoronoiClassifier3D(const int *regions)
void Initialize(vtkVoronoiClassifier3D *c)
bool IsInsideRegion(vtkIdType ptId)
bool IsSameRegion(vtkIdType ptId, vtkIdType neiId)
const vtkVoronoiSpoke * AddAdjacencyInformation(vtkVoronoiHull &hull, vtkVoronoiWheelsType &wheels, vtkVoronoiSpokesType &spokes, int &numSpokes, int &maxPoints, int &maxFaces)
ProduceWheelsAndSpokes(vtkVoronoiCore3D< TCompositor, TClassifier > *vc)
static void Execute(vtkVoronoiCore3D< TCompositor, TClassifier > *vc)
vtkVoronoiCore3D< TCompositor, TClassifier > * VC
void operator()(vtkIdType threadId, vtkIdType endThreadId)
void operator()(vtkIdType threadId, vtkIdType endThreadId)
vtkVoronoiCore3D< TCompositor, TClassifier > * VC
static std::unique_ptr< TopologicalMerge > Execute(vtkVoronoiCore3D< TCompositor, TClassifier > *vc)
vtkIdType GetNumberOfMergedPoints()
Methods related to merging coincident points.
TopologicalMerge(vtkVoronoiCore3D< TCompositor, TClassifier > *vc)
Typedefs and classes in support of the adjacency graph.
int vtkIdType
Definition vtkType.h:354
std::vector< vtkVoronoi2DLocalData< TCompositor, TClassifier > * > ThreadMapType
The thread map keeps track of the thread local data across all computing threads.
std::vector< vtkVoronoi3DLocalData< TCompositor, TClassifier > * > ThreadMapType
The thread map keeps track of the thread local data across all computing threads.
std::vector< vtkIdType > vtkBatchIdsType
Keep track of batches of generating points.
std::vector< vtkIdType > vtkMergeMapType
When merging points, the merge map is a vector that maps global tile/hull vertex ids (which contain d...
std::vector< vtkIdType > vtkVoronoiWheelsType
std::vector< vtkVoronoiSpoke > vtkVoronoiSpokesType
The vtkVoronoiWheelsType vector is used to keep track of the number of spokes (and equivalently,...
std::vector< vtkVoronoiMergeTuple3D > vtkMergeTuples3DType
#define vtkAlgorithm