VTK  9.6.20260214
vtkVoronoiTile Class Reference

The convex polygon tile class proper. More...

#include <vtkVoronoiTile.h>

Collaboration diagram for vtkVoronoiTile:
[legend]

Public Member Functions

 vtkVoronoiTile ()
 Constructor.
 
void Initialize (vtkIdType genPtId, const double genPt[3], double bds[4])
 Method to initiate the construction of the polygon.
 
void Initialize (vtkIdType ptId, const double x[3], vtkPoints *pts, vtkIdType nPts, const vtkIdType *p)
 Initialize with a convex polygon.
 
ClipIntersectionStatus Clip (vtkIdType neiPtId, const double neiPt[2])
 Method to clip the current convex tile with a plane defined by a neighboring point.
 
double GetCircumFlower2 ()
 Methods to determine whether a point x[2] is within the Voronoi flower, or Voronoi circumflower.
 
bool InCircumFlower (double r2)
 
bool InFlower (const double x[2])
 
void UpdatePetals (double cf2)
 
vtkDoubleArrayGetPetals ()
 
void ProducePolyData (vtkPolyData *pd, vtkSpheres *spheres)
 Produce a vtkPolyData (and optional implicit function) from the current polygon.
 
void ProducePolyData (vtkPolyData *pd)
 
vtkIdType GetGeneratorPointId ()
 Obtain information about the generated tile.
 
double * GetGeneratorPosition ()
 
vtkIdType GetNumberOfPoints ()
 
const PointRingTypeGetPoints ()
 

Public Attributes

vtkIdType PtId
 
double X [3]
 
vtkIdType NumClips
 
double PruneTolerance
 
PointRingType Points
 
PointRingType NewPoints
 

Protected Member Functions

PointRingIterator Next (PointRingIterator &itr)
 
ClipIntersectionStatus IntersectWithLine (double origin[2], double normal[2], vtkIdType neiPtId)
 
void ComputeCircumFlower ()
 

Protected Attributes

double Tol
 
double Tol2
 
bool RecomputeCircumFlower
 
bool RecomputePetals
 
double CircumFlower2
 
double MinRadius2
 
double MaxRadius2
 
std::vector< vtkTilePoint * > SortP
 
vtkSmartPointer< vtkDoubleArrayPetals
 

Detailed Description

The convex polygon tile class proper.

provide core 2D Voronoi tile generation capabilities

Since it is a supporting class, it is lightweight and not a subclass of vtkObject.

This lightweight, supporting class is used to generate a convex polygon (or tile) from repeated half-space clipping operations (i.e., generate a Voronoi tile). It also keeps track of the Voronoi flower and circumflower (i.e., the radius of security). These are used to determine whether a clipping operation will intersect the current Voronoi polygon.

The algorithm proceeds as follows. A generating point is placed within an initial, convex bounding box (i.e., this is the starting Voronoi tile). The hull is then repeatedly clipped by lines positioned at the halfway points between neighboring points, with each line's normal pointing in the direction of the edge connecting the generator point to a neighboring point.

The Voronoi tile class is represented by a counterclockwise ordered list of points. This also implicitly defines the Voronoi tile edges that form the polygon. In addition, the neighboring point ids–those which generated each polygon edge–are also maintained. This neighboring point information enables the production of topological constructs such as the Voronoi adjacency graph, which supports topological analysis capabilities.

Tolerancing capabilities are built into this class. The relative "PruneTolerance" is used to discard clipping nicks–that is, clipping lines that barely intersect (i.e., graze) the tile. By pruning (or discarding) small hull facets, the numerical stability of the tile generation process is significantly improved. Note that the PruneTolerance is relative, it is multiplied by a representative length of the tile; therefore it is adaptive to tile size.

Note
The tile is constructed in the x-y plane.
See also
vtkVoronoiHull vtkVoronoiCore2D vtkVoronoiFlower2D vtkVoronoiCore3D vtkVoronoi3D vtkGeneralizedSurfaceNets3D
Tests:
vtkVoronoiTile (Tests)

Definition at line 95 of file vtkVoronoiTile.h.

Constructor & Destructor Documentation

◆ vtkVoronoiTile()

vtkVoronoiTile::vtkVoronoiTile ( )
inline

Constructor.

After instantiation, for each point, make sure to initialize the tile with Initialize().

Definition at line 102 of file vtkVoronoiTile.h.

Member Function Documentation

◆ Initialize() [1/2]

void vtkVoronoiTile::Initialize ( vtkIdType genPtId,
const double genPt[3],
double bds[4] )

Method to initiate the construction of the polygon.

Define the generator point id and position, and an initial bounding box in which to place the generator.

◆ Initialize() [2/2]

void vtkVoronoiTile::Initialize ( vtkIdType ptId,
const double x[3],
vtkPoints * pts,
vtkIdType nPts,
const vtkIdType * p )

Initialize with a convex polygon.

The points must be in counterclockwise order (normal in the z-direction). Points must not be coincident. The polygon must be convex.

◆ Clip()

ClipIntersectionStatus vtkVoronoiTile::Clip ( vtkIdType neiPtId,
const double neiPt[2] )

Method to clip the current convex tile with a plane defined by a neighboring point.

The neighbor id and its position must not be coincident with the current generator point. This method does not take into account the Voronoi circumflower and flower.

◆ GetCircumFlower2()

double vtkVoronoiTile::GetCircumFlower2 ( )
inline

Methods to determine whether a point x[2] is within the Voronoi flower, or Voronoi circumflower.

(The Voronoi flower is the union of all Delaunay circles located at the tile points. The Voronoi circumflower is the 2*radius of the largest Delaunay circle.) These methods can be used to cull points which do not intersect the tile.

Definition at line 154 of file vtkVoronoiTile.h.

◆ InCircumFlower()

bool vtkVoronoiTile::InCircumFlower ( double r2)
inline

Definition at line 162 of file vtkVoronoiTile.h.

◆ InFlower()

bool vtkVoronoiTile::InFlower ( const double x[2])
inline

Definition at line 268 of file vtkVoronoiTile.h.

◆ UpdatePetals()

void vtkVoronoiTile::UpdatePetals ( double cf2)

◆ GetPetals()

vtkDoubleArray * vtkVoronoiTile::GetPetals ( )
inline

Definition at line 174 of file vtkVoronoiTile.h.

◆ ProducePolyData() [1/2]

void vtkVoronoiTile::ProducePolyData ( vtkPolyData * pd,
vtkSpheres * spheres )

Produce a vtkPolyData (and optional implicit function) from the current polygon.

If the spheres pointer is nullptr, it will not be generated. This method is typically used for debugging purposes.

◆ ProducePolyData() [2/2]

void vtkVoronoiTile::ProducePolyData ( vtkPolyData * pd)
inline

Definition at line 189 of file vtkVoronoiTile.h.

◆ GetGeneratorPointId()

vtkIdType vtkVoronoiTile::GetGeneratorPointId ( )
inline

Obtain information about the generated tile.

Note that in 2D, the number of points equals the number of convex polygon tile edges.

Definition at line 195 of file vtkVoronoiTile.h.

◆ GetGeneratorPosition()

double * vtkVoronoiTile::GetGeneratorPosition ( )
inline

Definition at line 196 of file vtkVoronoiTile.h.

◆ GetNumberOfPoints()

vtkIdType vtkVoronoiTile::GetNumberOfPoints ( )
inline

Definition at line 197 of file vtkVoronoiTile.h.

◆ GetPoints()

const PointRingType & vtkVoronoiTile::GetPoints ( )
inline

Definition at line 198 of file vtkVoronoiTile.h.

◆ Next()

PointRingIterator vtkVoronoiTile::Next ( PointRingIterator & itr)
inlineprotected

Definition at line 218 of file vtkVoronoiTile.h.

◆ IntersectWithLine()

ClipIntersectionStatus vtkVoronoiTile::IntersectWithLine ( double origin[2],
double normal[2],
vtkIdType neiPtId )
protected

◆ ComputeCircumFlower()

void vtkVoronoiTile::ComputeCircumFlower ( )
inlineprotected

Definition at line 248 of file vtkVoronoiTile.h.

Member Data Documentation

◆ PtId

vtkIdType vtkVoronoiTile::PtId

Definition at line 206 of file vtkVoronoiTile.h.

◆ X

double vtkVoronoiTile::X[3]

Definition at line 207 of file vtkVoronoiTile.h.

◆ NumClips

vtkIdType vtkVoronoiTile::NumClips

Definition at line 208 of file vtkVoronoiTile.h.

◆ PruneTolerance

double vtkVoronoiTile::PruneTolerance

Definition at line 209 of file vtkVoronoiTile.h.

◆ Points

PointRingType vtkVoronoiTile::Points

Definition at line 212 of file vtkVoronoiTile.h.

◆ NewPoints

PointRingType vtkVoronoiTile::NewPoints

Definition at line 213 of file vtkVoronoiTile.h.

◆ Tol

double vtkVoronoiTile::Tol
protected

Definition at line 231 of file vtkVoronoiTile.h.

◆ Tol2

double vtkVoronoiTile::Tol2
protected

Definition at line 232 of file vtkVoronoiTile.h.

◆ RecomputeCircumFlower

bool vtkVoronoiTile::RecomputeCircumFlower
protected

Definition at line 237 of file vtkVoronoiTile.h.

◆ RecomputePetals

bool vtkVoronoiTile::RecomputePetals
protected

Definition at line 238 of file vtkVoronoiTile.h.

◆ CircumFlower2

double vtkVoronoiTile::CircumFlower2
protected

Definition at line 239 of file vtkVoronoiTile.h.

◆ MinRadius2

double vtkVoronoiTile::MinRadius2
protected

Definition at line 240 of file vtkVoronoiTile.h.

◆ MaxRadius2

double vtkVoronoiTile::MaxRadius2
protected

Definition at line 241 of file vtkVoronoiTile.h.

◆ SortP

std::vector<vtkTilePoint*> vtkVoronoiTile::SortP
protected

Definition at line 242 of file vtkVoronoiTile.h.

◆ Petals

vtkSmartPointer<vtkDoubleArray> vtkVoronoiTile::Petals
protected

Definition at line 243 of file vtkVoronoiTile.h.


The documentation for this class was generated from the following file: