VTK  9.6.20260310
vtkVoronoiFlower3D.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
184
185#ifndef vtkVoronoiFlower3D_h
186#define vtkVoronoiFlower3D_h
187
188#include "vtkDataSetAlgorithm.h"
189#include "vtkFiltersMeshingModule.h" // For export macro
190#include "vtkIdTypeArray.h" // for PointsOfInterest
191#include "vtkSmartPointer.h" // For self-destructing data members
192#include "vtkStaticPointLocator.h" // For point locator
193#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
194
195VTK_ABI_NAMESPACE_BEGIN
196
197class VTKFILTERSMESHING_EXPORT VTK_MARSHALAUTO vtkVoronoiFlower3D : public vtkDataSetAlgorithm
198{
199public:
201
206 void PrintSelf(ostream& os, vtkIndent indent) override;
208
210
218 {
219 VORONOI = 0, // 3D Voronoi tessellation; output cells are polyhedra (vtkPolyhedron)
220 DELAUNAY = 1, // 3D Delaunay tesselation; cells are tetraheda (vtkTetra)
221 ADJACENCY_GRAPH = 2, // the graph edges connecting neighboring Voronoi hulls
223 3, // all polygonal faces including interior. Duplicate faces are not produced.
224 BOUNDARY = 4, // produce polygonal faces on the boundary of the Voronoi tessellation
226 5, // faces forming the surface net (i.e., faces on the boundaries between regions)
227 SPEED_TEST = 6, // no output, just compute Voronoi hulls (for performance testing)
228 };
229
237 vtkSetClampMacro(OutputType, int, VORONOI, SPEED_TEST);
238 vtkGetMacro(OutputType, int);
247
249
256 vtkSetClampMacro(Padding, double, 0.0001, 0.25);
257 vtkGetMacro(Padding, double);
259
261
270 vtkSetMacro(PassPointData, vtkTypeBool);
271 vtkGetMacro(PassPointData, vtkTypeBool);
272 vtkBooleanMacro(PassPointData, vtkTypeBool);
274
276
291 {
292 NO_CELL_SCALARS = 0, // Don't produce any cell scalars
293 POINT_IDS = 1, // Output cell scalars are the generating point id (default)
294 REGION_IDS = 2, // The region the cell primitives originated from (if region ids available)
295 NUMBER_FACES = 3, // The number of faces in the Voronoi hull
296 PRIM_IDS = 4, // The ids of the hull face primitives
298 5, // Scalars are the thread id used to produce output. This may change between runs.
299 RANDOM = 6 // Scalars are pseudo random numbers between [0,64).
300 };
301
302
304
313 vtkSetMacro(GenerateCellScalars, int);
314 vtkGetMacro(GenerateCellScalars, int);
323
325
332 vtkGetMacro(MergePoints, vtkTypeBool);
333 vtkSetMacro(MergePoints, vtkTypeBool);
334 vtkBooleanMacro(MergePoints, vtkTypeBool);
336
338
348 vtkIdType FindHull(double x[3]);
350
352
356 vtkSetClampMacro(PruneTolerance, double, 0.0, 0.5);
357 vtkGetMacro(PruneTolerance, double);
359
361
366 vtkStaticPointLocator* GetLocator() { return this->Locator; }
368
370
377 vtkSetMacro(Validate, vtkTypeBool);
378 vtkGetMacro(Validate, vtkTypeBool);
379 vtkBooleanMacro(Validate, vtkTypeBool);
381
383
391 vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
392 vtkGetMacro(BatchSize, unsigned int);
394
396
400 vtkGetMacro(BoundaryCapping, vtkTypeBool);
401 vtkSetMacro(BoundaryCapping, vtkTypeBool);
402 vtkBooleanMacro(BoundaryCapping, vtkTypeBool);
404
406
421 vtkSetClampMacro(PointOfInterest, vtkIdType, -1, VTK_ID_MAX);
422 vtkGetMacro(PointOfInterest, vtkIdType);
423 vtkSetObjectMacro(PointsOfInterest, vtkIdTypeArray);
424 vtkGetObjectMacro(PointsOfInterest, vtkIdTypeArray);
425 vtkSetClampMacro(MaximumNumberOfHullClips, vtkIdType, 1, VTK_ID_MAX);
426 vtkGetMacro(MaximumNumberOfHullClips, vtkIdType);
428
430
437 {
438 NO_POINT_SCALARS = 0, // Don't produce any point scalars (default)
439 FLOWER_RADII = 1 // Output point scalars are the Voronoi flower radii
440 };
441 vtkGetMacro(GeneratePointScalars, int);
443
449 int GetMaximumNumberOfPoints() { return this->MaximumNumberOfPoints; }
450
456 int GetMaximumNumberOfFaces() { return this->MaximumNumberOfFaces; }
457
463 int GetNumberOfThreads() { return this->NumberOfThreads; }
464
470 int GetNumberOfPrunes() { return this->NumberOfPrunes; }
471
476
481 template <typename T>
482 void UpdateExecutionInformation(T* voro);
483
484protected:
486 ~vtkVoronoiFlower3D() override = default;
487
488 // Satisfy pipeline-related API
490 int FillInputPortInformation(int port, vtkInformation* info) override;
491 int FillOutputPortInformation(int port, vtkInformation* info) override;
492
493private:
494 vtkVoronoiFlower3D(const vtkVoronoiFlower3D&) = delete;
495 void operator=(const vtkVoronoiFlower3D&) = delete;
496
497 int OutputType; // specification of the filter output
498 double Padding; // amount to pad out input points bounding box
499 vtkTypeBool Validate; // Choose to validate and repair output
500 vtkTypeBool PassPointData; // indicate whether to pass input point data to output
501 int GeneratePointScalars; // indicate whether point scalars are to be produced
502 int GenerateCellScalars; // indicate whether cell scalars are to be produced
503 vtkTypeBool MergePoints; // merge near coincident points or not
504 vtkIdType PointOfInterest; // specify a single input point to process
505 vtkSmartPointer<vtkIdTypeArray> PointsOfInterest; // list of points of interest
506 vtkIdType MaximumNumberOfHullClips; // limit the number of hull clips
507 vtkSmartPointer<vtkStaticPointLocator> Locator; // locator for finding proximal points
508 double PruneTolerance; // the prune spokes tolerance
509 unsigned int BatchSize; // process data in batches of specified size
510 vtkTypeBool BoundaryCapping; // cap the domain boundary if OutputType is SURFACE_NET
511
515 int NumberOfThreads; // report on the number of threads used during processing
516 int MaximumNumberOfPoints; // maximum number of points found in any hull
517 int MaximumNumberOfFaces; // maximum number of faces found in any hull
518 int NumberOfPrunes; // If spoke pruning is enabled, report number of pruning operations
519};
520
521//------------------------------------------------------------------------------
522template <typename T>
524{
525 this->NumberOfThreads = voro->GetNumberOfThreads();
526 this->MaximumNumberOfPoints = voro->GetMaximumNumberOfPoints();
527 this->MaximumNumberOfFaces = voro->GetMaximumNumberOfFaces();
528 this->NumberOfPrunes = voro->GetNumberOfPrunes();
529}
530
531VTK_ABI_NAMESPACE_END
532#endif
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Hold a reference to a vtkObjectBase instance.
quickly locate points in 3-space
void SetOutputTypeToBoundary()
Used to control the filter output.
int GetNumberOfPrunes()
Return the number of prunes performed during execution.
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
void SetGenerateCellScalarsToNumberFaces()
Indicate whether to create a cell scalar array as part of the output.
GenerateCellScalarsStrategy
Specify how to generate cell scalars for the outputs.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
~vtkVoronoiFlower3D() override=default
virtual void SetGenerateCellScalars(int)
Indicate whether to create a cell scalar array as part of the output.
void SetGenerateCellScalarsToPointIds()
Indicate whether to create a cell scalar array as part of the output.
vtkIdType FindHull(double x[3])
The following method–FindHull()–can be used to locate/query the Voronoi hull containing a point x (i....
OutputTypeOptions
Used to control the filter output.
GeneratePointScalarsStrategy
Used internally to generate point scalars for the output.
int GetNumberOfThreads()
Return the number of threads actually used during execution.
void SetOutputTypeToPolygonalComplex()
Used to control the filter output.
int GetMaximumNumberOfPoints()
Return the maximum number of points in any Voronoi hull.
void SetOutputTypeToVoronoi()
Used to control the filter output.
vtkStaticPointLocator * GetLocator()
Retrieve the internal locator to manually configure it, for example specifying the number of points p...
void SetGenerateCellScalarsToThreadIds()
Indicate whether to create a cell scalar array as part of the output.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
vtkMTimeType GetMTime() override
Get the MTime of this object also considering the locator.
virtual void SetOutputType(int)
Specify the type of output the filter creates.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called within ProcessRequest when a request asks the algorithm to do its work.
void UpdateExecutionInformation(T *voro)
Method used to update this filter's execution parameters after the internal, templated instance of vt...
void SetOutputTypeToDelaunay()
Used to control the filter output.
void SetGenerateCellScalarsToPrimIds()
Indicate whether to create a cell scalar array as part of the output.
void SetOutputTypeToAdjacencyGraph()
Used to control the filter output.
static vtkVoronoiFlower3D * New()
Standard methods for instantiation, type information, and printing.
int GetMaximumNumberOfFaces()
Return the maximum number of faces in any Voronoi hull.
void SetGenerateCellScalarsToRegionIds()
Indicate whether to create a cell scalar array as part of the output.
void SetGenerateCellScalarsToNone()
Indicate whether to create a cell scalar array as part of the output.
void SetOutputTypeToSurfaceNet()
Used to control the filter output.
void SetGenerateCellScalarsToRandom()
Indicate whether to create a cell scalar array as part of the output.
void SetOutputTypeToSpeedTest()
Used to control the filter output.
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:363
#define VTK_ID_MAX
Definition vtkType.h:367
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:318
#define VTK_INT_MAX
Definition vtkType.h:192
#define VTK_MARSHALAUTO