VTK  9.6.20260214
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
183
184#ifndef vtkVoronoiFlower3D_h
185#define vtkVoronoiFlower3D_h
186
187#include "vtkDataSetAlgorithm.h"
188#include "vtkFiltersMeshingModule.h" // For export macro
189#include "vtkIdTypeArray.h" // for PointsOfInterest
190#include "vtkSmartPointer.h" // For self-destructing data members
191#include "vtkStaticPointLocator.h" // For point locator
192#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
193
194VTK_ABI_NAMESPACE_BEGIN
195
196class VTKFILTERSMESHING_EXPORT VTK_MARSHALAUTO vtkVoronoiFlower3D : public vtkDataSetAlgorithm
197{
198public:
200
205 void PrintSelf(ostream& os, vtkIndent indent) override;
207
209
217 {
218 VORONOI = 0, // 3D Voronoi tessellation; output cells are polyhedra (vtkPolyhedron)
219 DELAUNAY = 1, // 3D Delaunay tesselation; cells are tetraheda (vtkTetra)
220 ADJACENCY_GRAPH = 2, // the graph edges connecting neighboring Voronoi hulls
222 3, // all polygonal faces including interior. Duplicate faces are not produced.
223 BOUNDARY = 4, // produce polygonal faces on the boundary of the Voronoi tessellation
225 5, // faces forming the surface net (i.e., faces on the boundaries between regions)
226 SPEED_TEST = 6, // no output, just compute Voronoi hulls (for performance testing)
227 };
228
236 vtkSetClampMacro(OutputType, int, VORONOI, SPEED_TEST);
237 vtkGetMacro(OutputType, int);
246
248
255 vtkSetClampMacro(Padding, double, 0.0001, 0.25);
256 vtkGetMacro(Padding, double);
258
260
269 vtkSetMacro(PassPointData, vtkTypeBool);
270 vtkGetMacro(PassPointData, vtkTypeBool);
271 vtkBooleanMacro(PassPointData, vtkTypeBool);
273
275
290 {
291 NO_CELL_SCALARS = 0, // Don't produce any cell scalars
292 POINT_IDS = 1, // Output cell scalars are the generating point id (default)
293 REGION_IDS = 2, // The region the cell primitives originated from (if region ids available)
294 NUMBER_FACES = 3, // The number of faces in the Voronoi hull
295 PRIM_IDS = 4, // The ids of the hull face primitives
297 5, // Scalars are the thread id used to produce output. This may change between runs.
298 RANDOM = 6 // Scalars are pseudo random numbers between [0,64).
299 };
300
301
303
312 vtkSetMacro(GenerateCellScalars, int);
313 vtkGetMacro(GenerateCellScalars, int);
322
324
331 vtkGetMacro(MergePoints, vtkTypeBool);
332 vtkSetMacro(MergePoints, vtkTypeBool);
333 vtkBooleanMacro(MergePoints, vtkTypeBool);
335
337
347 vtkIdType FindHull(double x[3]);
349
351
355 vtkSetClampMacro(PruneTolerance, double, 0.0, 0.5);
356 vtkGetMacro(PruneTolerance, double);
358
360
365 vtkStaticPointLocator* GetLocator() { return this->Locator; }
367
369
376 vtkSetMacro(Validate, vtkTypeBool);
377 vtkGetMacro(Validate, vtkTypeBool);
378 vtkBooleanMacro(Validate, vtkTypeBool);
380
382
390 vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
391 vtkGetMacro(BatchSize, unsigned int);
393
395
399 vtkGetMacro(BoundaryCapping, vtkTypeBool);
400 vtkSetMacro(BoundaryCapping, vtkTypeBool);
401 vtkBooleanMacro(BoundaryCapping, vtkTypeBool);
403
405
420 vtkSetClampMacro(PointOfInterest, vtkIdType, -1, VTK_ID_MAX);
421 vtkGetMacro(PointOfInterest, vtkIdType);
422 vtkSetObjectMacro(PointsOfInterest, vtkIdTypeArray);
423 vtkGetObjectMacro(PointsOfInterest, vtkIdTypeArray);
424 vtkSetClampMacro(MaximumNumberOfHullClips, vtkIdType, 1, VTK_ID_MAX);
425 vtkGetMacro(MaximumNumberOfHullClips, vtkIdType);
427
429
436 {
437 NO_POINT_SCALARS = 0, // Don't produce any point scalars (default)
438 FLOWER_RADII = 1 // Output point scalars are the Voronoi flower radii
439 };
440 vtkGetMacro(GeneratePointScalars, int);
442
448 int GetMaximumNumberOfPoints() { return this->MaximumNumberOfPoints; }
449
455 int GetMaximumNumberOfFaces() { return this->MaximumNumberOfFaces; }
456
462 int GetNumberOfThreads() { return this->NumberOfThreads; }
463
469 int GetNumberOfPrunes() { return this->NumberOfPrunes; }
470
475
480 template <typename T>
481 void UpdateExecutionInformation(T* voro);
482
483protected:
485 ~vtkVoronoiFlower3D() override = default;
486
487 // Satisfy pipeline-related API
489 int FillInputPortInformation(int port, vtkInformation* info) override;
490 int FillOutputPortInformation(int port, vtkInformation* info) override;
491
492private:
493 vtkVoronoiFlower3D(const vtkVoronoiFlower3D&) = delete;
494 void operator=(const vtkVoronoiFlower3D&) = delete;
495
496 int OutputType; // specification of the filter output
497 double Padding; // amount to pad out input points bounding box
498 vtkTypeBool Validate; // Choose to validate and repair output
499 vtkTypeBool PassPointData; // indicate whether to pass input point data to output
500 int GeneratePointScalars; // indicate whether point scalars are to be produced
501 int GenerateCellScalars; // indicate whether cell scalars are to be produced
502 vtkTypeBool MergePoints; // merge near coincident points or not
503 vtkIdType PointOfInterest; // specify a single input point to process
504 vtkSmartPointer<vtkIdTypeArray> PointsOfInterest; // list of points of interest
505 vtkIdType MaximumNumberOfHullClips; // limit the number of hull clips
506 vtkSmartPointer<vtkStaticPointLocator> Locator; // locator for finding proximal points
507 double PruneTolerance; // the prune spokes tolerance
508 unsigned int BatchSize; // process data in batches of specified size
509 vtkTypeBool BoundaryCapping; // cap the domain boundary if OutputType is SURFACE_NET
510
514 int NumberOfThreads; // report on the number of threads used during processing
515 int MaximumNumberOfPoints; // maximum number of points found in any hull
516 int MaximumNumberOfFaces; // maximum number of faces found in any hull
517 int NumberOfPrunes; // If spoke pruning is enabled, report number of pruning operations
518};
519
520//------------------------------------------------------------------------------
521template <typename T>
523{
524 this->NumberOfThreads = voro->GetNumberOfThreads();
525 this->MaximumNumberOfPoints = voro->GetMaximumNumberOfPoints();
526 this->MaximumNumberOfFaces = voro->GetMaximumNumberOfFaces();
527 this->NumberOfPrunes = voro->GetNumberOfPrunes();
528}
529
530VTK_ABI_NAMESPACE_END
531#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