VTK  9.6.20260312
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
197
198#ifndef vtkVoronoiFlower2D_h
199#define vtkVoronoiFlower2D_h
200
201#include "vtkAbstractTransform.h" // For transforming input points
202#include "vtkFiltersMeshingModule.h" // For export macro
203#include "vtkIdTypeArray.h" // for PointsOfInterest
204#include "vtkPolyDataAlgorithm.h"
205#include "vtkSmartPointer.h" // For self-destructing data members
206#include "vtkSpheres.h" // For Voronoi Flower
207#include "vtkStaticPointLocator2D.h" // For point locator
208#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
209
210VTK_ABI_NAMESPACE_BEGIN
213class vtkPointSet;
214class vtkSpheres;
215
216class VTKFILTERSMESHING_EXPORT VTK_MARSHALAUTO vtkVoronoiFlower2D : public vtkPolyDataAlgorithm
217{
218public:
220
225 void PrintSelf(ostream& os, vtkIndent indent) override;
227
239
241
263 vtkSetMacro(OutputType, int);
264 vtkGetMacro(OutputType, int);
271
273
284 vtkSetClampMacro(Padding, double, 0.0001, 0.25);
285 vtkGetMacro(Padding, double);
287
289
295 vtkSetMacro(PassPointData, vtkTypeBool);
296 vtkGetMacro(PassPointData, vtkTypeBool);
297 vtkBooleanMacro(PassPointData, vtkTypeBool);
299
301
308 {
309 NO_POINT_SCALARS = 0, // Don't produce any point scalars (default)
310 FLOWER_RADII = 1 // Output cell scalars are the generating point id
311 };
312 vtkGetMacro(GeneratePointScalars, int);
314
327 {
328 NO_CELL_SCALARS = 0, // Don't produce any cell scalars
329 POINT_IDS = 1, // Output cell scalars are the generating point id / tile id (default)
330 REGION_IDS = 2, // The region the cell primitives originated from (if region ids available)
331 NUMBER_SIDES = 3, // The number of edges in the Voronoi tile
332 PRIM_IDS = 4, // The ids of derivative primitives (e.g., Delaunay triangles)
334 5, // Scalars are the thread id used to produce output. This may change between runs.
335 RANDOM = 6 // Scalars are pseudo random numbers between [0,64).
336 };
337
339
348 vtkSetClampMacro(GenerateCellScalars, int, NO_CELL_SCALARS, RANDOM);
349 vtkGetMacro(GenerateCellScalars, int);
358
360
367 vtkGetMacro(MergePoints, vtkTypeBool);
368 vtkSetMacro(MergePoints, vtkTypeBool);
369 vtkBooleanMacro(MergePoints, vtkTypeBool);
371
373
386 vtkIdType FindTile(double x[3]);
387 void GetTileData(vtkIdType tileId, vtkPolyData* tileData);
389
391
400 vtkSetClampMacro(PruneTolerance, double, 0.0, 0.5);
401 vtkGetMacro(PruneTolerance, double);
403
405
412 vtkSetMacro(Validate, vtkTypeBool);
413 vtkGetMacro(Validate, vtkTypeBool);
414 vtkBooleanMacro(Validate, vtkTypeBool);
416
418
427 vtkSetSmartPointerMacro(Transform, vtkAbstractTransform);
428 vtkGetSmartPointerMacro(Transform, vtkAbstractTransform);
430
437
439
447 vtkSetClampMacro(ProjectionPlaneMode, int, XY_PLANE, BEST_FITTING_PLANE);
448 vtkGetMacro(ProjectionPlaneMode, int);
458
459
461
465 vtkGetMacro(BoundaryCapping, vtkTypeBool);
466 vtkSetMacro(BoundaryCapping, vtkTypeBool);
467 vtkBooleanMacro(BoundaryCapping, vtkTypeBool);
469
471
486 vtkSetClampMacro(PointOfInterest, vtkIdType, -1, VTK_ID_MAX);
487 vtkGetMacro(PointOfInterest, vtkIdType);
488 vtkSetObjectMacro(PointsOfInterest, vtkIdTypeArray);
489 vtkGetObjectMacro(PointsOfInterest, vtkIdTypeArray);
490 vtkSetClampMacro(MaximumNumberOfTileClips, vtkIdType, 1, VTK_ID_MAX);
491 vtkGetMacro(MaximumNumberOfTileClips, vtkIdType);
493
495
500 vtkStaticPointLocator2D* GetLocator() { return this->Locator; }
502
504
513 vtkSetMacro(GenerateVoronoiFlower, vtkTypeBool);
514 vtkGetMacro(GenerateVoronoiFlower, vtkTypeBool);
515 vtkBooleanMacro(GenerateVoronoiFlower, vtkTypeBool);
517
519
526 vtkGetObjectMacro(Spheres, vtkSpheres);
528
530
538 vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
539 vtkGetMacro(BatchSize, unsigned int);
541
543
548 int GetMaximumNumberOfPoints() { return this->MaximumNumberOfPoints; }
549 int GetMaximumNumberOfSides() { return this->MaximumNumberOfPoints; }
550 int GetMaximumNumberOfEdges() { return this->MaximumNumberOfPoints; }
552
558 int GetNumberOfThreads() { return this->NumberOfThreads; }
559
564
569 template <typename T>
570 void UpdateExecutionInformation(T* voro);
571
572protected:
574 ~vtkVoronoiFlower2D() override = default;
575
576 // Satisfy pipeline-related API
578 // Specify that the input is of type vtkPointSet
580
581private:
582 vtkVoronoiFlower2D(const vtkVoronoiFlower2D&) = delete;
583 void operator=(const vtkVoronoiFlower2D&) = delete;
584
585 int OutputType;
586 vtkTypeBool Validate;
587 double Padding;
588 vtkTypeBool PassPointData;
589 int GeneratePointScalars;
590 int GenerateCellScalars;
591 vtkTypeBool MergePoints; // merge near coincident points or not
592 int ProjectionPlaneMode; // selects the plane in 3D where the tessellation will be computed
595 vtkIdType PointOfInterest;
596 vtkSmartPointer<vtkIdTypeArray> PointsOfInterest; // list of points of interest
597 vtkIdType MaximumNumberOfTileClips;
598 vtkTypeBool GenerateVoronoiFlower;
600 double PruneTolerance;
601 unsigned int BatchSize;
602 vtkTypeBool BoundaryCapping; // cap the domain boundary if OutputType is SURFACE_NET
603
607 int NumberOfThreads; // report on the number of threads used during processing
608 int MaximumNumberOfPoints; // maximum number of points found in any hull
609 int NumberOfPrunes; // If spoke prining is enabled, report number of pruning operations
610};
611
612//------------------------------------------------------------------------------
613template <typename T>
615{
616 this->NumberOfThreads = voro->GetNumberOfThreads();
617 this->MaximumNumberOfPoints = voro->GetMaximumNumberOfPoints();
618 this->NumberOfPrunes = voro->GetNumberOfPrunes();
619}
620
621VTK_ABI_NAMESPACE_END
622#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 the type of output to generate: Voronoi, Delaunay, Voronoi+Delaunay, a surface net,...
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 the type of output to generate: Voronoi, Delaunay, Voronoi+Delaunay, a surface net,...
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 the type of output to generate: Voronoi, Delaunay, Voronoi+Delaunay, a surface net,...
void SetOutputTypeToDelaunay()
Control the type of output to generate: Voronoi, Delaunay, Voronoi+Delaunay, a surface net,...
void SetOutputTypeToSurfaceNet()
Control the type of output to generate: Voronoi, Delaunay, Voronoi+Delaunay, a surface net,...
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 the type of output to generate: Voronoi, Delaunay, Voronoi+Delaunay, a surface net,...
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