VTK
vtkVoronoi2D.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkVoronoi2D.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
94 #ifndef vtkVoronoi2D_h
95 #define vtkVoronoi2D_h
96 
97 #include "vtkFiltersCoreModule.h" // For export macro
98 #include "vtkPolyDataAlgorithm.h"
99 
102 class vtkPointSet;
103 class vtkSpheres;
104 
105 class VTKFILTERSCORE_EXPORT vtkVoronoi2D : public vtkPolyDataAlgorithm
106 {
107 public:
109 
112  static vtkVoronoi2D *New();
114  void PrintSelf(ostream& os, vtkIndent indent) override;
116 
118 
124  vtkSetClampMacro(Padding,double,0.001,0.25);
125  vtkGetMacro(Padding,double);
127 
129  {
130  NONE=0,
131  POINT_IDS=1,
132  THREAD_IDS=2
133  };
134 
136 
141  vtkSetMacro(GenerateScalars,int);
142  vtkGetMacro(GenerateScalars,int);
144  { this->SetGenerateScalars(NONE); }
146  { this->SetGenerateScalars(POINT_IDS); }
148  { this->SetGenerateScalars(THREAD_IDS); }
150 
152 
161  virtual void SetTransform(vtkAbstractTransform*);
162  vtkGetObjectMacro(Transform, vtkAbstractTransform);
164 
166  {
167  XY_PLANE=0,
168  SPECIFIED_TRANSFORM_PLANE=1,
169  BEST_FITTING_PLANE=2
170  };
171 
173 
181  vtkSetClampMacro(ProjectionPlaneMode,int,XY_PLANE,BEST_FITTING_PLANE);
182  vtkGetMacro(ProjectionPlaneMode,int);
184  { this->SetProjectionPlaneMode(XY_PLANE); }
186  { this->SetProjectionPlaneMode(SPECIFIED_TRANSFORM_PLANE); }
188  { this->SetProjectionPlaneMode(BEST_FITTING_PLANE); }
190 
192 
205  vtkSetClampMacro(PointOfInterest,vtkIdType,-1,VTK_ID_MAX);
206  vtkGetMacro(PointOfInterest,vtkIdType);
207  vtkSetClampMacro(MaximumNumberOfTileClips,vtkIdType,1,VTK_ID_MAX);
208  vtkGetMacro(MaximumNumberOfTileClips,vtkIdType);
210 
212 
218  { return this->Locator; }
220 
222 
231  vtkSetMacro(GenerateVoronoiFlower,vtkTypeBool);
232  vtkGetMacro(GenerateVoronoiFlower,vtkTypeBool);
233  vtkBooleanMacro(GenerateVoronoiFlower,vtkTypeBool);
235 
237 
244  vtkGetObjectMacro(Spheres,vtkSpheres);
246 
252  {return this->NumberOfThreadsUsed;}
253 
257  vtkMTimeType GetMTime() override;
258 
259 protected:
260  vtkVoronoi2D();
261  ~vtkVoronoi2D() override;
262 
264  double Padding;
265  double Tolerance;
266  int ProjectionPlaneMode; //selects the plane in 3D where the tessellation will be computed
274 
275  // Satisfy pipeline-related API
277  int FillInputPortInformation(int, vtkInformation*) override;
278 
279 private:
280  vtkVoronoi2D(const vtkVoronoi2D&) = delete;
281  void operator=(const vtkVoronoi2D&) = delete;
282 };
283 
284 #endif
vtkSpheres * Spheres
Definition: vtkVoronoi2D.h:273
int ProjectionPlaneMode
Definition: vtkVoronoi2D.h:266
implicit function for a set of spheres
Definition: vtkSpheres.h:42
Store vtkAlgorithm input/output information.
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:302
void SetProjectionPlaneModeToSpecifiedTransformPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
Definition: vtkVoronoi2D.h:185
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
vtkAbstractTransform * Transform
Definition: vtkVoronoi2D.h:268
vtkStaticPointLocator2D * Locator
Definition: vtkVoronoi2D.h:267
abstract class for specifying dataset behavior
Definition: vtkPointSet.h:42
int vtkIdType
Definition: vtkType.h:347
static vtkPolyDataAlgorithm * New()
int vtkTypeBool
Definition: vtkABI.h:69
vtkIdType MaximumNumberOfTileClips
Definition: vtkVoronoi2D.h:270
Superclass for algorithms that produce only polydata as output.
a simple class to control print indentation
Definition: vtkIndent.h:39
double Tolerance
Definition: vtkVoronoi2D.h:265
superclass for all geometric transformations
vtkTypeBool GenerateVoronoiFlower
Definition: vtkVoronoi2D.h:271
virtual vtkMTimeType GetMTime()
Return this object's modified time.
create 2D Voronoi convex tiling of input points
Definition: vtkVoronoi2D.h:105
int GetNumberOfThreadsUsed()
Return the number of threads actually used during execution.
Definition: vtkVoronoi2D.h:251
double Padding
Definition: vtkVoronoi2D.h:264
void SetProjectionPlaneModeToBestFittingPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
Definition: vtkVoronoi2D.h:187
#define VTK_ID_MAX
Definition: vtkType.h:351
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
Store zero or more vtkInformation instances.
vtkStaticPointLocator2D * GetLocator()
Retrieve the internal locator to manually configure it, for example specifying the number of points p...
Definition: vtkVoronoi2D.h:217
int NumberOfThreadsUsed
Definition: vtkVoronoi2D.h:272
void SetGenerateScalarsToNone()
Indicate whether to create a scalar array as part of the output.
Definition: vtkVoronoi2D.h:143
void SetProjectionPlaneModeToXYPlane()
Define the method to project the input 3D points into a 2D plane for tessellation.
Definition: vtkVoronoi2D.h:183
void SetGenerateScalarsToThreadIds()
Indicate whether to create a scalar array as part of the output.
Definition: vtkVoronoi2D.h:147
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkIdType PointOfInterest
Definition: vtkVoronoi2D.h:269
void SetGenerateScalarsToPointIds()
Indicate whether to create a scalar array as part of the output.
Definition: vtkVoronoi2D.h:145
quickly locate points in 2-space
Transform
Definition: ADIOSDefs.h:40