VTK  9.6.20260213
vtkVoronoiCore3D< TCompositor, TClassifier > Class Template Reference

========================================================================= The templated, core Voronoi class. More...

#include <vtkVoronoiCore3D.h>

Inheritance diagram for vtkVoronoiCore3D< TCompositor, TClassifier >:
[legend]
Collaboration diagram for vtkVoronoiCore3D< TCompositor, TClassifier >:
[legend]

Classes

struct  ProduceWheelsAndSpokes
 Produce the global adjacency graph / wheels and spokes data structure. More...
 
struct  TopologicalMerge
 Functor class used to topologically merge (nearly) coincident points. More...
 

Public Member Functions

vtkVoronoiAdjacencyGraphGetAdjacencyGraph ()
 Obtain the adjacency graph (wheel & spokes data structure).
 
 vtkVoronoiCore3D (vtkAlgorithm *filter, vtkVoronoiBatchManager &batcher, vtkStaticPointLocator *loc, vtkPoints *inPts, double padding, vtkIdType maxClips, bool validate, double pruneTol, TCompositor *comp, TClassifier *c)
 Constructor.
 
vtkIdType GetNumberOfPoints ()
 Convenience methods to retrieve the number of input points, and the raw double* points array.
 
const double * GetPoints ()
 
int GetNumberOfThreads ()
 Access the local thread data produced by execution of the filter.
 
vtkVoronoi3DLocalData< TCompositor, TClassifier > * GetThreadData (int threadNum)
 Access the local thread data produced by execution of the filter.
 
int GetMaximumNumberOfPoints ()
 Obtain information about the execution of the Voronoi algorithm.
 
int GetMaximumNumberOfFaces ()
 Obtain information about the execution of the Voronoi algorithm.
 
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.
 
void Reduce ()
 Core vtkSMPTools methods to process hulls in parallel, and produce thread local output data.
 

Static Public Member Functions

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.
 

Public Attributes

TCompositor Compositor
 The compositor enables this vtkVoronoiCore3D templated class to be used in different applications.
 
TClassifier Classifier
 This templated class is used to extend the API of this vtkVoronoiCore3D class, to implement the spoke classification process, to clone copies in multiple threads, and to initialize the classification instances.
 
vtkVoronoiBatchManager Batcher
 Controls processing of batches of generating points.
 
ThreadMapType< TCompositor, TClassifier > ThreadMap
 
vtkVoronoiAdjacencyGraph Graph
 This is used to create the spokes and wheels adjacency graph used to validate the tessellation and produce a Delaunay triangulation.
 
vtkAlgorithmFilter
 Used for controlling filter abort and accessing filter information.
 

Detailed Description

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
class vtkVoronoiCore3D< TCompositor, TClassifier >

========================================================================= The templated, core Voronoi class.

provide core 3D Voronoi tessellation capabilities

It is a lightweight supporting class (i.e., not a subclass of vtkObject) meant to be used by specialized algorithms requiring Voronoi and/or Delaunay capabilities.

Note: the template argument TCompositor is used to control what information is extracted during tessellation. Different using filters will define and extract information relevant to their application needs. This is accomplished by defining different compositing classes. TClassifier is used to classify the spokes connecting neighborhood points, which to the the dual property, classifies the hull faces.

The Voronoi tessellation is a common computational tool used in a variety of applications ranging from tessellating points, mesh generation, surface reconstruction, materials anaysis, and contouring (surface nets). It can also be the basis for computing its dual construct, the Delaunay triangulation, also used in wide-ranging applications with significant impacts. This templated class provides core 3D Voronoi tessellation capabilities, including implementation of fast parallel algorithms, which can be used by other classes to create specialized Voronoi-based algorithms.

The Voronoi tessellation is a tessellation of space, where each Voronoi n-dimensional tile (in 3D, a polyhedron, also referred to as a point hull, or more simply, hull) represents the region nearest to one of the input points. Under the hood, the vtkVoronoiCore3D class depends on the vtkVoronoiHull class to construct hulls. The vtkVoronoiCore3D provides a framework for parallel, SMP shared memory algorithms, which can be specialized to meet the needs of algorithms requiring Voronoi and/or Delaunay-based capabilities. Specialization of the algorithm is via a TCompositor template argument, which controls what information is extracted from each hull as it is processed, and how the information is combined (during a compositing process) to produce output.

Publications are in preparation to describe the algorithm. Conceptually, the algorithm is meshless, meaning that each input point and its associated hull is processed independently. However, methods are provided to transform the output into a fully-connected, valid, conformal mesh if required. To summarize: in parallel, each (generating) input point is associated with an initial Voronoi hull, which is simply the bounding box of the input point set. A spatial locator is then used to identify nearby points: each neighbor in turn generates a clipping plane positioned halfway between the generating point and the neighboring point, and orthogonal to the line connecting them. Clips are performed by evaluationg the vertices of the convex Voronoi hull as being on either side (inside,outside) of the clip plane. If an intersection with the Voronoi hull is found, the portion of the hull "outside" the clip plane is discarded, producing a new convex face, resulting in a new convex, Voronoi hull. As each clip occurs, the Voronoi "Flower" error metric (the union of Delunay spheres) is compared to the extent of the region containing the neighboring clip points. The clip region (along with the points contained in it) is grown by careful expansion, When the Voronoi circumflower is contained within the clip region, the algorithm terminates and the Voronoi hull is output. Once complete, it is possible to construct the Delaunay triangulation from the Voronoi tessellation. Note that topological and geometric information can be used to generate a valid triangulation (e.g., merging coincident points and validating topology).

This class can also construct a Voronoi adjacency graph, composed of edges (the spokes) that connect each Voronoi hull generating point (the wheels) with their face neighbors. The adjacency graph is a powerful data representation that captures proximal neighborhood information. It has many practical applications such as shortest path computation.

An implementation note: this class is implemented using a parallel algorithm. The output is invariant no matter what order the the threads execute, i.e., the construction of geometric primitives (Voronoi cells, adjacency graphs, etc.) is identical no matter the number of threads used, or execution order. Also, note the correspondance between input generating point of ptId and hullId, ptId produces hulls of hullId, where ptId == hullId. This means for debugging purposes, picking output primitives with POINT_IDS enabled provides a means to select the original generating hull.

This class is templated, enabling specialized capabilities depending on the using algorithm / filter. Using templating, it is possible to extract different information from each Voronoi hull as it is constructed, compositing this information later to produce different types of output.

Warning
Coincident input points are discarded. The Voronoi tessellation requires unique input points.
This approach implements an embarrassingly parallel, meshless algorithm. At the core of the algorithm a locator is used to determine points close to a specified position. Currently, a vtkStaticPointLocator is used because it is both threaded (when constructed) and supports thread-safe queries. While other locators could be used in principal, they must support thread-safe operations. This is done by defining and implementing an iterator that enables traversal of neighborhood points surrounding each generating point.
This class has been threaded with vtkSMPTools. Using TBB or other non-sequential type (set in the CMake variable VTK_SMP_IMPLEMENTATION_TYPE) may improve performance significantly. For example, in recent benchmarks using a 48-thread desktop, one million points can be processed in well under one second elapsed time.
See also
vtkVoronoiHull vtkVoronoi3D vtkGeneralizedSurfaceNets3D
Tests:
vtkVoronoiCore3D (Tests)

Definition at line 242 of file vtkVoronoiCore3D.h.

Constructor & Destructor Documentation

◆ vtkVoronoiCore3D()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
vtkVoronoiCore3D< TCompositor, TClassifier >::vtkVoronoiCore3D ( vtkAlgorithm * filter,
vtkVoronoiBatchManager & batcher,
vtkStaticPointLocator * loc,
vtkPoints * inPts,
double padding,
vtkIdType maxClips,
bool validate,
double pruneTol,
TCompositor * comp,
TClassifier * c )

Constructor.

This is typically not directly invoked by the user. The Execute() method is preferred.

Member Function Documentation

◆ Execute()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
static std::unique_ptr< vtkVoronoiCore3D< TCompositor, TClassifier > > 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 )
static

A factory method to conveniently instantiate and execute the algorithm.

This class should be executed using this Execute() method. It returns a unique pointer to an instance of the VoronoiCore3D algorithm. After execution, methods on the instance can be invoked to retrieve relevant information. Note also that a Voronoi compositor should also be provided to this Execute() method. It will contain output as well. Input to the method includes the input (double) points, a prebuilt static point locator, the initial hull bounds padding, a limit on the number of clips each hull can perform, and an optional VTK filter (for controlling execution abort). An optional (non-null) compositor and/or classifier can be provided which is used to initialize the compositors and/or classifiers in the various threads (for example, providing region ids). Finally, methods to control degenerate faces (i.e., validation / spoke pruning and tolerance) are provided.

◆ GetNumberOfThreads()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
int vtkVoronoiCore3D< TCompositor, TClassifier >::GetNumberOfThreads ( )
inline

Access the local thread data produced by execution of the filter.

This includes the compositing data. The data is only available after the Execute() method has been invoked.

Definition at line 271 of file vtkVoronoiCore3D.h.

◆ GetThreadData()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
vtkVoronoi3DLocalData< TCompositor, TClassifier > * vtkVoronoiCore3D< TCompositor, TClassifier >::GetThreadData ( int threadNum)
inline

Access the local thread data produced by execution of the filter.

This includes the compositing data. The data is only available after the Execute() method has been invoked.

Definition at line 272 of file vtkVoronoiCore3D.h.

◆ GetMaximumNumberOfPoints()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
int vtkVoronoiCore3D< TCompositor, TClassifier >::GetMaximumNumberOfPoints ( )
inline

Obtain information about the execution of the Voronoi algorithm.

This includes the maximum number of faces found in any hull; the maximum number of points found in any hull; and the number of prunes performed to remove degeneracies.

Definition at line 285 of file vtkVoronoiCore3D.h.

◆ GetMaximumNumberOfFaces()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
int vtkVoronoiCore3D< TCompositor, TClassifier >::GetMaximumNumberOfFaces ( )
inline

Obtain information about the execution of the Voronoi algorithm.

This includes the maximum number of faces found in any hull; the maximum number of points found in any hull; and the number of prunes performed to remove degeneracies.

Definition at line 286 of file vtkVoronoiCore3D.h.

◆ GetNumberOfPrunes()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
int vtkVoronoiCore3D< TCompositor, TClassifier >::GetNumberOfPrunes ( )
inline

Obtain information about the execution of the Voronoi algorithm.

This includes the maximum number of faces found in any hull; the maximum number of points found in any hull; and the number of prunes performed to remove degeneracies.

Definition at line 287 of file vtkVoronoiCore3D.h.

◆ GetAdjacencyGraph()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
vtkVoronoiAdjacencyGraph & vtkVoronoiCore3D< TCompositor, TClassifier >::GetAdjacencyGraph ( )
inline

Obtain the adjacency graph (wheel & spokes data structure).

This is constructed during algorithm execution.

Definition at line 294 of file vtkVoronoiCore3D.h.

◆ Initialize()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
void vtkVoronoiCore3D< TCompositor, TClassifier >::Initialize ( )

Core vtkSMPTools methods to process hulls in parallel, and produce thread local output data.

Other compositing data is captured as well. Note that the threading occurs over batches of points. (Note: this is left public for access by vtkSMPTools. These are generally not invoked directly by the user.)

◆ operator()()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
void vtkVoronoiCore3D< TCompositor, TClassifier >::operator() ( vtkIdType batchId,
vtkIdType endBatchId )

Core vtkSMPTools methods to process hulls in parallel, and produce thread local output data.

Other compositing data is captured as well. Note that the threading occurs over batches of points. (Note: this is left public for access by vtkSMPTools. These are generally not invoked directly by the user.)

◆ Reduce()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
void vtkVoronoiCore3D< TCompositor, TClassifier >::Reduce ( )

Core vtkSMPTools methods to process hulls in parallel, and produce thread local output data.

Other compositing data is captured as well. Note that the threading occurs over batches of points. (Note: this is left public for access by vtkSMPTools. These are generally not invoked directly by the user.)

◆ GetNumberOfPoints()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
vtkIdType vtkVoronoiCore3D< TCompositor, TClassifier >::GetNumberOfPoints ( )
inline

Convenience methods to retrieve the number of input points, and the raw double* points array.

Invoke this only after execution.

Definition at line 351 of file vtkVoronoiCore3D.h.

◆ GetPoints()

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
const double * vtkVoronoiCore3D< TCompositor, TClassifier >::GetPoints ( )
inline

Definition at line 352 of file vtkVoronoiCore3D.h.

Member Data Documentation

◆ Compositor

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
TCompositor vtkVoronoiCore3D< TCompositor, TClassifier >::Compositor

The compositor enables this vtkVoronoiCore3D templated class to be used in different applications.

It supports parallel gather/compute of specified information on a hull-by-hull basis, which can then combined/composited to produce output. Users of this class must define their own compositor.

Definition at line 324 of file vtkVoronoiCore3D.h.

◆ Classifier

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
TClassifier vtkVoronoiCore3D< TCompositor, TClassifier >::Classifier

This templated class is used to extend the API of this vtkVoronoiCore3D class, to implement the spoke classification process, to clone copies in multiple threads, and to initialize the classification instances.

Definition at line 331 of file vtkVoronoiCore3D.h.

◆ Batcher

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
vtkVoronoiBatchManager vtkVoronoiCore3D< TCompositor, TClassifier >::Batcher

Controls processing of batches of generating points.

Thread local data is is available after generating the hulls.

Definition at line 337 of file vtkVoronoiCore3D.h.

◆ ThreadMap

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
ThreadMapType<TCompositor, TClassifier> vtkVoronoiCore3D< TCompositor, TClassifier >::ThreadMap

Definition at line 338 of file vtkVoronoiCore3D.h.

◆ Graph

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
vtkVoronoiAdjacencyGraph vtkVoronoiCore3D< TCompositor, TClassifier >::Graph

This is used to create the spokes and wheels adjacency graph used to validate the tessellation and produce a Delaunay triangulation.

Note that if an "empty" classifier is used, the adjacency graph is empty.

Definition at line 345 of file vtkVoronoiCore3D.h.

◆ Filter

template<class TCompositor, class TClassifier = vtkVoronoiClassifier3D>
vtkAlgorithm* vtkVoronoiCore3D< TCompositor, TClassifier >::Filter

Used for controlling filter abort and accessing filter information.

If nullptr, then filter abort checking is disabled.

Definition at line 358 of file vtkVoronoiCore3D.h.


The documentation for this class was generated from the following file: