VTK  9.3.20240419
vtkVoronoi2D.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
82 #ifndef vtkVoronoi2D_h
83 #define vtkVoronoi2D_h
84 
85 #include "vtkFiltersCoreModule.h" // For export macro
86 #include "vtkPolyDataAlgorithm.h"
87 
88 VTK_ABI_NAMESPACE_BEGIN
91 class vtkPointSet;
92 class vtkSpheres;
93 
94 class VTKFILTERSCORE_EXPORT vtkVoronoi2D : public vtkPolyDataAlgorithm
95 {
96 public:
98 
101  static vtkVoronoi2D* New();
103  void PrintSelf(ostream& os, vtkIndent indent) override;
105 
107 
113  vtkSetClampMacro(Padding, double, 0.001, 0.25);
114  vtkGetMacro(Padding, double);
116 
118  {
119  NONE = 0,
120  POINT_IDS = 1,
121  THREAD_IDS = 2
122  };
123 
125 
130  vtkSetMacro(GenerateScalars, int);
131  vtkGetMacro(GenerateScalars, int);
132  void SetGenerateScalarsToNone() { this->SetGenerateScalars(NONE); }
133  void SetGenerateScalarsToPointIds() { this->SetGenerateScalars(POINT_IDS); }
134  void SetGenerateScalarsToThreadIds() { this->SetGenerateScalars(THREAD_IDS); }
136 
138 
148  vtkGetObjectMacro(Transform, vtkAbstractTransform);
150 
152  {
153  XY_PLANE = 0,
154  SPECIFIED_TRANSFORM_PLANE = 1,
155  BEST_FITTING_PLANE = 2
156  };
157 
159 
167  vtkSetClampMacro(ProjectionPlaneMode, int, XY_PLANE, BEST_FITTING_PLANE);
168  vtkGetMacro(ProjectionPlaneMode, int);
169  void SetProjectionPlaneModeToXYPlane() { this->SetProjectionPlaneMode(XY_PLANE); }
171  {
172  this->SetProjectionPlaneMode(SPECIFIED_TRANSFORM_PLANE);
173  }
175  {
176  this->SetProjectionPlaneMode(BEST_FITTING_PLANE);
177  }
179 
181 
194  vtkSetClampMacro(PointOfInterest, vtkIdType, -1, VTK_ID_MAX);
195  vtkGetMacro(PointOfInterest, vtkIdType);
196  vtkSetClampMacro(MaximumNumberOfTileClips, vtkIdType, 1, VTK_ID_MAX);
197  vtkGetMacro(MaximumNumberOfTileClips, vtkIdType);
199 
201 
206  vtkStaticPointLocator2D* GetLocator() { return this->Locator; }
208 
210 
219  vtkSetMacro(GenerateVoronoiFlower, vtkTypeBool);
220  vtkGetMacro(GenerateVoronoiFlower, vtkTypeBool);
221  vtkBooleanMacro(GenerateVoronoiFlower, vtkTypeBool);
223 
225 
232  vtkGetObjectMacro(Spheres, vtkSpheres);
234 
239  int GetNumberOfThreadsUsed() { return this->NumberOfThreadsUsed; }
240 
244  vtkMTimeType GetMTime() override;
245 
246 protected:
248  ~vtkVoronoi2D() override;
249 
251  double Padding;
252  double Tolerance;
253  int ProjectionPlaneMode; // selects the plane in 3D where the tessellation will be computed
261 
262  // Satisfy pipeline-related API
265 
266 private:
267  vtkVoronoi2D(const vtkVoronoi2D&) = delete;
268  void operator=(const vtkVoronoi2D&) = delete;
269 };
270 
271 VTK_ABI_NAMESPACE_END
272 #endif
superclass for all geometric transformations
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
Superclass for algorithms that produce only polydata as output.
implicit function for a set of spheres
Definition: vtkSpheres.h:32
quickly locate points in 2-space
create 2D Voronoi convex tiling of input points
Definition: vtkVoronoi2D.h:95
vtkSpheres * Spheres
Definition: vtkVoronoi2D.h:260
vtkStaticPointLocator2D * GetLocator()
Retrieve the internal locator to manually configure it, for example specifying the number of points p...
Definition: vtkVoronoi2D.h:206
void SetGenerateScalarsToNone()
Indicate whether to create a scalar array as part of the output.
Definition: vtkVoronoi2D.h:132
vtkAbstractTransform * Transform
Definition: vtkVoronoi2D.h:255
double Tolerance
Definition: vtkVoronoi2D.h:252
void SetProjectionPlaneModeToSpecifiedTransformPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
Definition: vtkVoronoi2D.h:170
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
vtkIdType MaximumNumberOfTileClips
Definition: vtkVoronoi2D.h:257
~vtkVoronoi2D() override
void SetGenerateScalarsToPointIds()
Indicate whether to create a scalar array as part of the output.
Definition: vtkVoronoi2D.h:133
void SetProjectionPlaneModeToBestFittingPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
Definition: vtkVoronoi2D.h:174
void SetProjectionPlaneModeToXYPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
Definition: vtkVoronoi2D.h:169
int NumberOfThreadsUsed
Definition: vtkVoronoi2D.h:259
static vtkVoronoi2D * New()
Standard methods for instantiation, type information, and printing.
vtkStaticPointLocator2D * Locator
Definition: vtkVoronoi2D.h:254
int GetNumberOfThreadsUsed()
Return the number of threads actually used during execution.
Definition: vtkVoronoi2D.h:239
vtkIdType PointOfInterest
Definition: vtkVoronoi2D.h:256
int ProjectionPlaneMode
Definition: vtkVoronoi2D.h:253
virtual void SetTransform(vtkAbstractTransform *)
Set / get the transform which is applied to points to generate a 2D problem.
vtkMTimeType GetMTime() override
Get the MTime of this object also considering the locator.
int FillInputPortInformation(int, vtkInformation *) override
Fill the input port information objects for this algorithm.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkTypeBool GenerateVoronoiFlower
Definition: vtkVoronoi2D.h:258
void SetGenerateScalarsToThreadIds()
Indicate whether to create a scalar array as part of the output.
Definition: vtkVoronoi2D.h:134
double Padding
Definition: vtkVoronoi2D.h:251
@ Transform
Definition: vtkX3D.h:41
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315
#define VTK_ID_MAX
Definition: vtkType.h:319
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270