VTK  9.5.20250905
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
87#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
88
89VTK_ABI_NAMESPACE_BEGIN
92class vtkPointSet;
93class vtkSpheres;
94
95class VTKFILTERSCORE_EXPORT VTK_MARSHALAUTO vtkVoronoi2D : public vtkPolyDataAlgorithm
96{
97public:
99
102 static vtkVoronoi2D* New();
104 void PrintSelf(ostream& os, vtkIndent indent) override;
106
108
114 vtkSetClampMacro(Padding, double, 0.001, 0.25);
115 vtkGetMacro(Padding, double);
117
119 {
120 NONE = 0,
121 POINT_IDS = 1,
122 THREAD_IDS = 2
123 };
124
126
131 vtkSetMacro(GenerateScalars, int);
132 vtkGetMacro(GenerateScalars, int);
133 void SetGenerateScalarsToNone() { this->SetGenerateScalars(NONE); }
134 void SetGenerateScalarsToPointIds() { this->SetGenerateScalars(POINT_IDS); }
135 void SetGenerateScalarsToThreadIds() { this->SetGenerateScalars(THREAD_IDS); }
137
139
149 vtkGetObjectMacro(Transform, vtkAbstractTransform);
151
153 {
154 XY_PLANE = 0,
155 SPECIFIED_TRANSFORM_PLANE = 1,
156 BEST_FITTING_PLANE = 2
157 };
158
160
168 vtkSetClampMacro(ProjectionPlaneMode, int, XY_PLANE, BEST_FITTING_PLANE);
169 vtkGetMacro(ProjectionPlaneMode, int);
170 void SetProjectionPlaneModeToXYPlane() { this->SetProjectionPlaneMode(XY_PLANE); }
172 {
173 this->SetProjectionPlaneMode(SPECIFIED_TRANSFORM_PLANE);
174 }
176 {
177 this->SetProjectionPlaneMode(BEST_FITTING_PLANE);
178 }
180
182
195 vtkSetClampMacro(PointOfInterest, vtkIdType, -1, VTK_ID_MAX);
196 vtkGetMacro(PointOfInterest, vtkIdType);
197 vtkSetClampMacro(MaximumNumberOfTileClips, vtkIdType, 1, VTK_ID_MAX);
198 vtkGetMacro(MaximumNumberOfTileClips, vtkIdType);
200
202
207 vtkStaticPointLocator2D* GetLocator() { return this->Locator; }
209
211
220 vtkSetMacro(GenerateVoronoiFlower, vtkTypeBool);
221 vtkGetMacro(GenerateVoronoiFlower, vtkTypeBool);
222 vtkBooleanMacro(GenerateVoronoiFlower, vtkTypeBool);
224
226
233 vtkGetObjectMacro(Spheres, vtkSpheres);
235
240 int GetNumberOfThreadsUsed() { return this->NumberOfThreadsUsed; }
241
246
247protected:
249 ~vtkVoronoi2D() override;
250
252 double Padding;
253 double Tolerance;
254 int ProjectionPlaneMode; // selects the plane in 3D where the tessellation will be computed
262
263 // Satisfy pipeline-related API
266
267private:
268 vtkVoronoi2D(const vtkVoronoi2D&) = delete;
269 void operator=(const vtkVoronoi2D&) = delete;
270};
271
272VTK_ABI_NAMESPACE_END
273#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
vtkSpheres * Spheres
void SetGenerateScalarsToNone()
Indicate whether to create a 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...
vtkAbstractTransform * Transform
void SetProjectionPlaneModeToSpecifiedTransformPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
vtkIdType MaximumNumberOfTileClips
~vtkVoronoi2D() override
void SetGenerateScalarsToPointIds()
Indicate whether to create a scalar array as part of the output.
void SetProjectionPlaneModeToBestFittingPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
void SetProjectionPlaneModeToXYPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
int NumberOfThreadsUsed
static vtkVoronoi2D * New()
Standard methods for instantiation, type information, and printing.
vtkStaticPointLocator2D * Locator
int GetNumberOfThreadsUsed()
Return the number of threads actually used during execution.
vtkIdType PointOfInterest
int ProjectionPlaneMode
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
void SetGenerateScalarsToThreadIds()
Indicate whether to create a scalar array as part of the output.
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:332
#define VTK_ID_MAX
Definition vtkType.h:336
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287
#define VTK_MARSHALAUTO