VTK  9.6.20260215
vtkVoronoiFlower2D.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
185
186#ifndef vtkVoronoiFlower2D_h
187#define vtkVoronoiFlower2D_h
188
189#include "vtkAbstractTransform.h" // For transforming input points
190#include "vtkFiltersMeshingModule.h" // For export macro
191#include "vtkIdTypeArray.h" // for PointsOfInterest
192#include "vtkPolyDataAlgorithm.h"
193#include "vtkSmartPointer.h" // For self-destructing data members
194#include "vtkSpheres.h" // For Voronoi Flower
195#include "vtkStaticPointLocator2D.h" // For point locator
196#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
197
198VTK_ABI_NAMESPACE_BEGIN
201class vtkPointSet;
202class vtkSpheres;
203
204class VTKFILTERSMESHING_EXPORT VTK_MARSHALAUTO vtkVoronoiFlower2D : public vtkPolyDataAlgorithm
205{
206public:
208
213 void PrintSelf(ostream& os, vtkIndent indent) override;
215
226
228
246 vtkSetMacro(OutputType, int);
247 vtkGetMacro(OutputType, int);
253
255
266 vtkSetClampMacro(Padding, double, 0.0001, 0.25);
267 vtkGetMacro(Padding, double);
269
271
277 vtkSetMacro(PassPointData, vtkTypeBool);
278 vtkGetMacro(PassPointData, vtkTypeBool);
279 vtkBooleanMacro(PassPointData, vtkTypeBool);
281
283
290 {
291 NO_POINT_SCALARS = 0, // Don't produce any point scalars (default)
292 FLOWER_RADII = 1 // Output cell scalars are the generating point id
293 };
294 vtkGetMacro(GeneratePointScalars, int);
296
309 {
310 NO_CELL_SCALARS = 0, // Don't produce any cell scalars
311 POINT_IDS = 1, // Output cell scalars are the generating point id / tile id (default)
312 REGION_IDS = 2, // The region the cell primitives originated from (if region ids available)
313 NUMBER_SIDES = 3, // The number of edges in the Voronoi tile
314 PRIM_IDS = 4, // The ids of derivative primitives (e.g., Delaunay triangles)
316 5, // Scalars are the thread id used to produce output. This may change between runs.
317 RANDOM = 6 // Scalars are pseudo random numbers between [0,64).
318 };
319
321
330 vtkSetClampMacro(GenerateCellScalars, int, NO_CELL_SCALARS, RANDOM);
331 vtkGetMacro(GenerateCellScalars, int);
340
342
349 vtkGetMacro(MergePoints, vtkTypeBool);
350 vtkSetMacro(MergePoints, vtkTypeBool);
351 vtkBooleanMacro(MergePoints, vtkTypeBool);
353
355
368 vtkIdType FindTile(double x[3]);
369 void GetTileData(vtkIdType tileId, vtkPolyData* tileData);
371
373
382 vtkSetClampMacro(PruneTolerance, double, 0.0, 0.5);
383 vtkGetMacro(PruneTolerance, double);
385
387
394 vtkSetMacro(Validate, vtkTypeBool);
395 vtkGetMacro(Validate, vtkTypeBool);
396 vtkBooleanMacro(Validate, vtkTypeBool);
398
400
409 vtkSetSmartPointerMacro(Transform, vtkAbstractTransform);
410 vtkGetSmartPointerMacro(Transform, vtkAbstractTransform);
412
419
421
429 vtkSetClampMacro(ProjectionPlaneMode, int, XY_PLANE, BEST_FITTING_PLANE);
430 vtkGetMacro(ProjectionPlaneMode, int);
440
441
443
458 vtkSetClampMacro(PointOfInterest, vtkIdType, -1, VTK_ID_MAX);
459 vtkGetMacro(PointOfInterest, vtkIdType);
460 vtkSetObjectMacro(PointsOfInterest, vtkIdTypeArray);
461 vtkGetObjectMacro(PointsOfInterest, vtkIdTypeArray);
462 vtkSetClampMacro(MaximumNumberOfTileClips, vtkIdType, 1, VTK_ID_MAX);
463 vtkGetMacro(MaximumNumberOfTileClips, vtkIdType);
465
467
472 vtkStaticPointLocator2D* GetLocator() { return this->Locator; }
474
476
485 vtkSetMacro(GenerateVoronoiFlower, vtkTypeBool);
486 vtkGetMacro(GenerateVoronoiFlower, vtkTypeBool);
487 vtkBooleanMacro(GenerateVoronoiFlower, vtkTypeBool);
489
491
498 vtkGetObjectMacro(Spheres, vtkSpheres);
500
502
510 vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
511 vtkGetMacro(BatchSize, unsigned int);
513
515
520 int GetMaximumNumberOfPoints() { return this->MaximumNumberOfPoints; }
521 int GetMaximumNumberOfSides() { return this->MaximumNumberOfPoints; }
522 int GetMaximumNumberOfEdges() { return this->MaximumNumberOfPoints; }
524
530 int GetNumberOfThreads() { return this->NumberOfThreads; }
531
536
541 template <typename T>
542 void UpdateExecutionInformation(T* voro);
543
544protected:
546 ~vtkVoronoiFlower2D() override = default;
547
548 // Satisfy pipeline-related API
550 // Specify that the input is of type vtkPointSet
552
553private:
554 vtkVoronoiFlower2D(const vtkVoronoiFlower2D&) = delete;
555 void operator=(const vtkVoronoiFlower2D&) = delete;
556
557 int OutputType;
558 vtkTypeBool Validate;
559 double Padding;
560 vtkTypeBool PassPointData;
561 int GeneratePointScalars;
562 int GenerateCellScalars;
563 vtkTypeBool MergePoints; // merge near coincident points or not
564 int ProjectionPlaneMode; // selects the plane in 3D where the tessellation will be computed
567 vtkIdType PointOfInterest;
568 vtkSmartPointer<vtkIdTypeArray> PointsOfInterest; // list of points of interest
569 vtkIdType MaximumNumberOfTileClips;
570 vtkTypeBool GenerateVoronoiFlower;
572 double PruneTolerance;
573 unsigned int BatchSize;
574
578 int NumberOfThreads; // report on the number of threads used during processing
579 int MaximumNumberOfPoints; // maximum number of points found in any hull
580 int NumberOfPrunes; // If spoke prining is enabled, report number of pruning operations
581};
582
583//------------------------------------------------------------------------------
584template <typename T>
586{
587 this->NumberOfThreads = voro->GetNumberOfThreads();
588 this->MaximumNumberOfPoints = voro->GetMaximumNumberOfPoints();
589 this->NumberOfPrunes = voro->GetNumberOfPrunes();
590}
591
592VTK_ABI_NAMESPACE_END
593#endif
superclass for all geometric transformations
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.
concrete class for storing a set of points
Definition vtkPointSet.h:98
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
implicit function for a set of spheres
Definition vtkSpheres.h:32
quickly locate points in 2-space
void SetProjectionPlaneModeToSpecifiedTransformPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
~vtkVoronoiFlower2D() override=default
void SetGenerateCellScalarsToRandom()
Indicate whether to create a cell scalar array as part of the output.
void SetOutputTypeToVoronoiAndDelaunay()
Control whether to produce an output Voronoi tessellation and/or an output Delaunay triangulation.
void SetGenerateCellScalarsToPrimIds()
Indicate whether to create a cell scalar array as part of the output.
virtual void SetGenerateCellScalars(int)
Indicate whether to create a cell scalar array as part of the output.
OutputTypeOptions
Used to control filter output.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int GetMaximumNumberOfSides()
Return the maximum number of sides across all Voronoi tiles.
int GetMaximumNumberOfPoints()
Return the maximum number of sides across all Voronoi tiles.
int GetNumberOfThreads()
Return the number of threads actually used during execution.
void SetOutputTypeToSpeedTest()
Control whether to produce an output Voronoi tessellation and/or an output Delaunay triangulation.
int GetMaximumNumberOfEdges()
Return the maximum number of sides across all Voronoi tiles.
void SetGenerateCellScalarsToPointIds()
Indicate whether to create a cell scalar array as part of the output.
static vtkVoronoiFlower2D * New()
Standard methods for instantiation, type information, and printing.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
GeneratePointScalarsStrategy
Used internally to generate point scalars for the output.
virtual void SetOutputType(int)
Control whether to produce an output Voronoi tessellation and/or an output Delaunay triangulation.
void SetOutputTypeToDelaunay()
Control whether to produce an output Voronoi tessellation and/or an output Delaunay triangulation.
int FillInputPortInformation(int, vtkInformation *) override
Fill the input port information objects for this algorithm.
void SetProjectionPlaneModeToBestFittingPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
vtkIdType FindTile(double x[3])
The following methods - FindTile() and GetTileData() - can be used to locate/query the tile containin...
void GetTileData(vtkIdType tileId, vtkPolyData *tileData)
The following methods - FindTile() and GetTileData() - can be used to locate/query the tile containin...
vtkMTimeType GetMTime() override
Get the MTime of this object also considering the locator.
void SetGenerateCellScalarsToNumberOfSides()
Indicate whether to create a cell scalar array as part of the output.
void SetOutputTypeToVoronoi()
Control whether to produce an output Voronoi tessellation and/or an output Delaunay triangulation.
void UpdateExecutionInformation(T *voro)
Method used to update this filter's execution parameters after the internal, templated instance of vt...
void SetProjectionPlaneModeToXYPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
void SetGenerateCellScalarsToThreadIds()
Indicate whether to create a cell scalar array as part of the output.
void SetGenerateCellScalarsToRegionIds()
Indicate whether to create a cell scalar array as part of the output.
vtkStaticPointLocator2D * GetLocator()
Retrieve the internal locator to manually configure it, for example specifying the number of points p...
void SetGenerateCellScalarsToNone()
Indicate whether to create a cell scalar array as part of the output.
GenerateCellScalarsStrategy
Specify how to generate cell scalars for the outputs.
virtual void SetProjectionPlaneMode(int)
Define the method to project the input 3D points into a 2D plane for tessellation.
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