VTK  9.6.20260215
vtkVoronoiHull.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
45
46#ifndef vtkVoronoiHull_h
47#define vtkVoronoiHull_h
48
49#include "vtkDoubleArray.h" // Support double points
50#include "vtkFiltersMeshingModule.h" // For export macro
51#include "vtkMath.h" // Euclidean norm
52#include "vtkVoronoiCore.h" // Enum: clip intersection status
53
54#include <vector> // array support
55
56VTK_ABI_NAMESPACE_BEGIN
57
58class vtkPolyData;
59
68{
71};
72
73//======= Define the convex polyhedron class used to produce Voronoi hulls.
74
87struct VTKFILTERSMESHING_EXPORT vtkHullPoint
88{
89 double X[3]; // position
90 double Val; // evaluated value against a half-space clipping plane
91 double R2; // Voronoi petal radius
92 int PtMap; // renumber points that are actually used by valid faces
93 int Faces[3]; // the three faces defining this point
94 ProcessingStatus Status; // the status of the point
95
96 // methods to define new point
97 vtkHullPoint(double x, double y, double z)
98 : X{ x, y, z }
99 , Val(0)
100 , R2(0.0)
101 , PtMap(-1)
102 , Faces{ -1, -1, -1 }
104 {
105 }
106 vtkHullPoint(double x[3])
107 : X{ x[0], x[1], x[2] }
108 , Val(0)
109 , R2(0.0)
110 , PtMap(-1)
111 , Faces{ -1, -1, -1 }
113 {
114 }
115
116 // Replace a deleted point with a new point.
117 void Replace(double x, double y, double z)
118 {
119 this->X[0] = x;
120 this->X[1] = y;
121 this->X[2] = z;
122 this->Val = 0.0;
123 this->R2 = 0.0;
124 this->PtMap = (-1);
125 this->Faces[0] = (-1);
126 this->Faces[1] = (-1);
127 this->Faces[2] = (-1);
128 this->Status = ProcessingStatus::Valid;
129 }
130}; // vtkHullPoint
131
137struct VTKFILTERSMESHING_EXPORT vtkHullFace
138{
139 vtkIdType NeiId; // the neighboring generator point that produced this face
140 int NumPts; // the number of point ids defining this face
141 int Points; // the offset into FacePoints listing the point ids of this face
142 int AllocSize; // the number of slots allocated for the Points
143 int InsertPos; // the position to add a new point in the Points array
144 ProcessingStatus Status; // the status of this face
145
146 // Construct a polyhedron face
148 : NeiId(neiId)
149 , NumPts(0)
150 , Points(-1)
151 , AllocSize(0)
153 {
154 }
155
156 // Replace a deleted face with a new face. Allocated memory for point ids is
157 // untouched.
158 void Replace(vtkIdType neiPtId)
159 {
160 this->NeiId = neiPtId;
161 this->Status = ProcessingStatus::Valid;
162 }
163}; // vtkHullFace
164
182{
183 vtkIdType V0; // min edge vertex id
184 vtkIdType V1; // max edge vertex id
185 vtkIdType Id; // point id of inserted point
186 vtkIdType LoopIdx; // order of the point around the new face.
187 vtkIdType Faces[2]; // the two faces using the edge.
188
190 : V0(v0)
191 , V1(v1)
192 , Id(-1)
193 , LoopIdx(-1)
194 , Faces{ -1, -1 }
195 {
196 if (this->V0 > this->V1)
197 {
198 std::swap(this->V0, this->V1);
199 }
200 }
201 bool operator==(const vtkHullEdgeTuple& et) const
202 {
203 return this->V0 == et.V0 && this->V1 == et.V1;
204 }
205 bool operator!=(const vtkHullEdgeTuple& et) const
206 {
207 return this->V0 != et.V0 || this->V1 != et.V1;
208 }
209 bool IsEdge(vtkIdType v0, vtkIdType v1) const
210 {
211 if (v0 < v1) // ordered properly
212 {
213 return (this->V0 == v0 && this->V1 == v1);
214 }
215 else // swap comparison required
216 {
217 return (this->V0 == v1 && this->V1 == v0);
218 }
219 }
220}; // vtkHullEdgeTuple
221
226class VTKFILTERSMESHING_EXPORT vtkVoronoiHull
227{
228public:
230 std::vector<vtkHullPoint>; // the polyhedral points (and associated attributes)
231 using vtkHullFacesArray = std::vector<vtkHullFace>; // a list of polyhedral faces
232 using FacePointsArray = std::vector<int>; // the list of points (by id) defining the faces
233 using FaceScratchArray = std::vector<int>; // temporary face point ids array
234 using InsertedEdgePointsArray = std::vector<vtkHullEdgeTuple>; // collect edge intersection points
235
241 : PtId(-1)
242 , X{ 0, 0, 0 }
243 , NumClips(0)
244 , PruneTolerance(1.0e-13)
246 , RecomputePetals(true)
247 , CircumFlower2(0.0)
248 , MinRadius2(0.0)
249 , MaxRadius2(0.0)
250 {
251 // Preallocate some space
252 this->Points.reserve(256);
253 this->Faces.reserve(256);
254 this->FacePoints.reserve(2048);
255 this->InProcessPoints.reserve(256);
256 this->InProcessFaces.reserve(256);
257 this->DeletedPoints.Reserve(256);
258 this->DeletedFaces.Reserve(256);
259
260 // Supporting data structures
261 this->SortP.reserve(256);
263 this->Petals->SetNumberOfComponents(4); // x-y-z-R2
264 this->Petals->Allocate(256); // initial allocation
265 }
266
272 void Initialize(vtkIdType genPtId, const double genPt[3], double bds[6]);
273
281 ClipIntersectionStatus Clip(vtkIdType neiPtId, const double neiPt[3]);
282
290 double GetCircumFlower2() { return this->CircumFlower2; }
291 bool InCircumFlower(double r2) // radius**2 of point from generator
292 {
293 // Only recompute the circumflower if necessary; that is, when
294 // a maximal point is eliminated by a polyhedral plane clip.
295 if (this->RecomputeCircumFlower)
296 {
297 this->ComputeCircumFlower();
298 }
300 }
301 bool InFlower(const double x[3]);
302 void UpdatePetals(double cf2);
304 {
305 if (this->RecomputePetals)
306 {
307 this->UpdatePetals(this->CircumFlower2);
308 }
309 return (this->Petals->GetNumberOfTuples() > 0 ? this->Petals : nullptr);
310 }
311
316 void MapPoints();
317
323
329
330 // Data members purposely left public for using classes to extract
331 // information.
332
333 // Information used to define the polyhedron- its generating point id and
334 // position, plus region classification. Indicate whether degenerate faces
335 // (i.e., those having ~zero area) can be deleted (i.e., pruned).
336 vtkIdType PtId; // Generating point id
337 double X[3]; // Generating point position
338 vtkIdType NumClips; // The total number of clip operations since Initialize()
339 double PruneTolerance; // Specify the prune tolerance
340
341 // Support normal joggle in the case of degeneracies.
343 void BumpNormal(int bumpNum, double normal[3], double bumpNormal[3]);
344 void BumpOrigin(int bumpNum, double origin[3], double bumpOrigin[3]);
345
346 // These data members represent the constructed polyhedron.
347 vtkIdType NumPts; // The number of valid points in the points array
348 vtkHullPointsArray Points; // Array of points defining this polyhedron
349 vtkIdType NumFaces; // The number of valid faces in the faces array
350 vtkHullFacesArray Faces; // A list of faces forming this polyhedron
351
358 {
359 std::vector<int> Stack;
360 DeletionStack() { this->Stack.reserve(256); }
361 void Push(int id) { this->Stack.emplace_back(id); }
362 int Pop()
363 {
364 vtkIdType id = this->Stack.back();
365 this->Stack.pop_back();
366 return id;
367 }
368 bool IsEmpty() { return (this->Stack.empty()); }
369 void Clear() { this->Stack.clear(); }
370 void Reserve(int sze) { this->Stack.reserve(sze); }
371 }; // DeletionStack
372
376
381 {
382 std::function<void(vtkVoronoiHull& hull, int, int, int)> Function;
386
387 // Constructor
388 FaceOperation(int faceId)
389 : FaceId(faceId)
390 {
391 }
392
393 // Support std::find()
394 bool operator==(int faceId) const { return (this->FaceId == faceId); }
395 };
396
403 struct FaceProcessingArray : public std::vector<FaceOperation>
404 {
405 bool AddFace(vtkVoronoiHull* hull, int faceId)
406 {
407 if (std::find(this->begin(), this->end(), faceId) == this->end())
408 {
409 FaceOperation& faceOp = this->emplace_back(faceId);
410 int numInts = hull->EvaluateFace(hull->GetFace(faceId), faceOp.StartIdx, faceOp.NumKeptPts);
411 if (numInts <= 0)
412 {
414 }
415 else if (numInts == 2)
416 {
418 }
419 else
420 {
421 return false;
422 }
423 }
424 return true; // successfully added the face for processing
425 }
426 };
427
431 struct PointProcessingArray : public std::vector<int>
432 {
433 // Add a non-degenerate point, and connected faces, for
434 // processing.
435 bool AddPoint(vtkVoronoiHull* hull, vtkHullPoint& point, int ptId)
436 {
437 this->emplace_back(ptId);
438 if (hull->InProcessFaces.AddFace(hull, point.Faces[0]) &&
439 hull->InProcessFaces.AddFace(hull, point.Faces[1]) &&
440 hull->InProcessFaces.AddFace(hull, point.Faces[2]))
441 {
442 return true;
443 }
444 else
445 {
446 return false;
447 }
448 }
449 };
450
454 vtkIdType AddNewPoint(double x, double y, double z);
455 vtkIdType AddNewPoint(double x[3]) { return this->AddNewPoint(x[0], x[1], x[2]); }
456 void DeletePoint(vtkIdType ptId);
457 void SetPointFaces(vtkIdType pId, vtkIdType f0, vtkIdType f1, vtkIdType f2);
458
462 int AddNewFace(vtkIdType npts, vtkIdType neiPtId);
463 vtkHullFace* GetFace(int faceId) { return &(this->Faces[faceId]); }
469 {
470 this->FacePoints[face->Points + face->InsertPos++] = ptId;
471 }
472
476 void AddNthFacePoint(vtkHullFace* face, int idx, vtkIdType ptId)
477 {
478 this->FacePoints[face->Points + idx] = ptId;
479 }
480
483 int GetFacePoint(vtkHullFace* face, int ptNum) { return this->FacePoints[face->Points + ptNum]; }
484
488 void RebuildFacePoints(vtkHullFace* face, FaceScratchArray& idsBuffer);
489
504 int EvaluateFace(vtkHullFace* face, int& startIdx, int& numKeptPts);
505
513 void DeleteFace(int faceId, int vtkNotUsed(startIdx), int vtkNotUsed(numKeptPts));
514
522 void RebuildFace(int faceId, int startIdx, int numKeptPts);
523
529 int IntersectFaceEdge(int faceId, int p0, int p1);
530
535 void AllocatePointIds(int npts, vtkHullFace& face);
536
537protected:
541 FacePointsArray FacePoints; // Point ids used to define faces
542 int MIN_POINTIDS_ALLOC = 10; // Minimum buffer allocation for face point ids
543
544 // Keep track of deleted points and faces so their memory can be reused. This acts as a
545 // sort of poor person's memory pool.
548
549 // Used to process and track faces and points affected by a plane clip operation.
550 PointProcessingArray InProcessPoints; // Points affected by current clip operation
551 FaceProcessingArray InProcessFaces; // Faces affected by current clip operation
552 InsertedEdgePointsArray InsertedEdgePoints; // New points generated on intersected edges
553 FaceScratchArray FaceIdsBuffer; // An internal buffer used to rebuild faces after clipping
554
555 // Indicate whether the Voronoi circumflower needs recomputing, and
556 // keep track of the current circumflower and related information.
557 void ComputeCircumFlower();
563 std::vector<vtkHullPoint*> SortP; // Points sorted on radius**2
564 vtkSmartPointer<vtkDoubleArray> Petals; // Flower petals w/ radii > shell radius
565
566 // Internal methods for processing clip operations of various types.
571 void Clear()
572 {
573 this->NumPts = 0;
574 this->Points.clear();
575 this->NumFaces = 0;
576 this->Faces.clear();
577 this->FacePoints.clear();
578 this->DeletedPoints.Clear();
579 this->DeletedFaces.Clear();
580 }
581
586 ClipIntersectionStatus IntersectWithPlane(double origin[3], double normal[3], vtkIdType neiPtId);
587}; // vtkVoronoiHull
588
589// In the following, inlined methods for performance
590
591//------------------------------------------------------------------------------
593{
594 this->Points[pId].Faces[0] = f0;
595 this->Points[pId].Faces[1] = f1;
596 this->Points[pId].Faces[2] = f2;
597}
598
599//------------------------------------------------------------------------------
600inline vtkIdType vtkVoronoiHull::AddNewPoint(double x, double y, double z)
601{
602 vtkIdType id;
603 // If there are no empty slots in the points array, then allocate a new point.
604 if (this->DeletedPoints.IsEmpty())
605 {
606 id = this->Points.size();
607 this->Points.emplace_back(x, y, z);
608 this->Points.back().R2 = vtkMath::Distance2BetweenPoints(this->Points.back().X, this->X);
609 }
610 else // otherwise, replace a previously deleted point and reuse its memory.
611 {
612 id = this->DeletedPoints.Pop();
613 this->Points[id].Replace(x, y, z);
614 this->Points[id].R2 = vtkMath::Distance2BetweenPoints(this->Points[id].X, this->X);
615 }
616
617 this->NumPts++;
618 return id;
619}
620
621//------------------------------------------------------------------------------
623{
624 // The circumflower may need recomputing if this point is an extreme point.
625 if ((4 * this->Points[ptId].R2) >= this->CircumFlower2)
626 {
627 this->RecomputeCircumFlower = true;
628 }
629
630 this->DeletedPoints.Push(ptId);
631 this->Points[ptId].Status = ProcessingStatus::Deleted;
632 this->NumPts--;
633}
634
635//------------------------------------------------------------------------------
637{
638 vtkIdType id;
639 // If there are no empty slots in the faces array, than allocate a new face.
640 if (this->DeletedFaces.IsEmpty())
641 {
642 id = this->Faces.size();
643 this->Faces.emplace_back(neiPtId);
644 this->AllocatePointIds(npts, this->Faces.back());
645 }
646 else
647 {
648 id = this->DeletedFaces.Pop();
649 this->Faces[id].Replace(neiPtId);
650 this->AllocatePointIds(npts, this->Faces[id]);
651 }
652
653 this->NumFaces++;
654 return id;
655}
656
657//------------------------------------------------------------------------------
658inline int vtkVoronoiHull::EvaluateFace(vtkHullFace* face, int& startIdx, int& numKeptPts)
659{
660 startIdx = numKeptPts = 0;
661 int npts = face->NumPts, numEdgeInts = 0;
662 int ip, p0, p1;
663 double val0, val1;
664
665 for (int i = 0; i < npts; ++i)
666 {
667 p0 = this->GetFacePoint(face, i);
668 ip = ((i + 1) == npts ? 0 : (i + 1));
669 p1 = this->GetFacePoint(face, ip);
670
671 val0 = this->Points[p0].Val;
672 val1 = this->Points[p1].Val;
673
674 if (val0 <= 0.0)
675 {
676 numKeptPts++;
677 }
678
679 if (val0 > 0 && val1 <= 0)
680 {
681 numEdgeInts++;
682 startIdx = i;
683 }
684 else if (val0 <= 0 && val1 > 0)
685 {
686 numEdgeInts++;
687 }
688 }
689
690 return numEdgeInts;
691}
692
693//------------------------------------------------------------------------------
695{
696 face.NumPts = npts;
697 face.InsertPos = 0;
698
699 // See if allocation is necessary. Otherwise use previous.
700 if (npts > face.AllocSize)
701 {
702 int size = (npts > MIN_POINTIDS_ALLOC ? npts : MIN_POINTIDS_ALLOC);
703 vtkIdType offset = this->FacePoints.size();
704 this->FacePoints.insert(this->FacePoints.end(), size, (-1));
705 face.AllocSize = size;
706 face.Points = offset;
707 }
708}
709
710//------------------------------------------------------------------------------
712{
713 // Make sure space has been allocated for rebuilt point ids.
714 int npts = static_cast<int>(idsBuffer.size());
715 if (npts > face->AllocSize)
716 {
717 // Need to realloc space for point ids.
718 this->AllocatePointIds(npts, *face);
719 }
720
721 // Copy buffer into face point ids.
722 for (face->InsertPos = 0; face->InsertPos < npts; ++face->InsertPos)
723 {
724 this->FacePoints[face->Points + face->InsertPos] = idsBuffer[face->InsertPos];
725 }
726 face->NumPts = npts;
727}
728
729//------------------------------------------------------------------------------
730inline void vtkVoronoiHull::DeleteFace(int faceId, int, int)
731{
732 this->DeletedFaces.Push(faceId);
733 this->Faces[faceId].Status = ProcessingStatus::Deleted;
734 this->NumFaces--;
735}
736
737//------------------------------------------------------------------------------
739{
740 // Compute the circumflower, and compute some info about
741 // the flower radii.
744
745 // Determine the circumflower and minimal sphere radius by
746 // checking against each of the flower petals.
747 for (const auto& pt : this->Points)
748 {
749 if (pt.Status == ProcessingStatus::Valid)
750 {
751 this->MinRadius2 = std::min(pt.R2, this->MinRadius2);
752 this->MaxRadius2 = std::max(pt.R2, this->MaxRadius2);
753 }
754 }
755 this->CircumFlower2 = (4.0 * this->MaxRadius2); // (2*(max petal radius))**2
756 this->RecomputeCircumFlower = false; // circumflower is up to date
757}
758
759//------------------------------------------------------------------------------
760inline bool vtkVoronoiHull::InFlower(const double x[3])
761{
762 // Check against the flower petals
763 for (const auto& pt : this->Points)
764 {
765 if (pt.Status == ProcessingStatus::Valid)
766 {
767 double r2 = vtkMath::Distance2BetweenPoints(x, pt.X);
768 if (r2 <= pt.R2)
769 {
770 return true;
771 }
772 }
773 } // for all valid points
774 // Point not in the flower since it's not in any petals
775 return false;
776}
777
778//------------------------------------------------------------------------------
780{
781 // Renumber the output points. Note that associated Faces should use the
782 // PtMap id to ensure the the point connectivity ids are contiguous.
783 vtkIdType id = 0;
784 for (auto& pitr : this->Points)
785 {
786 if (pitr.Status == ProcessingStatus::Valid)
787 {
788 pitr.PtMap = id++;
789 }
790 }
791}
792
793//------------------------------------------------------------------------------
794// Update the flower petals which are passed off to the locator.
795// Only petals which exend past the minimal radius of the shell
796// request are added to the list of petals. It is presumed that
797// UpdateCircumFlower() has been invoked previously.
798inline void vtkVoronoiHull::UpdatePetals(double cf2)
799{
800 // If the radii of the flower spheres (petals) is highly variable (which
801 // occurs when the spacing of points is highly variable), then there is
802 // likely a lot of empty search space. Only add flower petals which extend
803 // past the outer shell request boundary. These petals are used to further
804 // limit the point search space.
805 this->Petals->Reset();
806 this->RecomputePetals = false; // petals will be updated in the following
807
808 constexpr double SphereRatio = 2.0, SphereRatio2 = SphereRatio * SphereRatio;
809 if (this->MinRadius2 > 0 && ((this->MaxRadius2 / this->MinRadius2) < SphereRatio2))
810 {
811 return; // it's not worth using the petals
812 }
813
814 // Empirically determined
815 constexpr double LargeSphereRatio = 0.25;
816 int maxLargeSpheres = LargeSphereRatio * this->NumPts;
817
818 this->SortP.clear();
819 double minR2 = VTK_FLOAT_MAX, maxR2 = VTK_FLOAT_MIN;
820 for (auto& pt : this->Points)
821 {
822 if (pt.Status == ProcessingStatus::Valid)
823 {
824 // (2*R)**2 >= shell request radius**2
825 if ((4 * pt.R2) >= cf2)
826 {
827 minR2 = std::min(pt.R2, minR2);
828 maxR2 = std::max(pt.R2, maxR2);
829 this->SortP.emplace_back(&pt);
830 }
831 }
832 }
833
834 if (static_cast<int>(this->SortP.size()) > maxLargeSpheres || (maxR2 / minR2) < SphereRatio2)
835 {
836 return; // it's not worth using the petals
837 }
838 else // sort from large spheres to small
839 {
840 std::sort(this->SortP.begin(), this->SortP.end(),
841 [](const vtkHullPoint* p0, const vtkHullPoint* p1) { return (p0->R2 > p1->R2); });
842 for (auto& pt : this->SortP)
843 {
844 this->Petals->InsertNextTuple4(pt->X[0], pt->X[1], pt->X[2], pt->R2);
845 }
846 } // it's worth using large sphere culling
847} // UpdatePetals
848
849VTK_ABI_NAMESPACE_END
850#endif
851// VTK-HeaderTest-Exclude: vtkVoronoiHull.h
RealT r2
Definition PyrC2Basis.h:20
dynamic, self-adjusting array of double
static ReturnTypeT Distance2BetweenPoints(const TupleRangeT1 &p1, const TupleRangeT2 &p2)
Compute distance squared between two points p1 and p2.
Definition vtkMath.h:2065
concrete dataset represents vertices, lines, polygons, and triangle strips
Hold a reference to a vtkObjectBase instance.
static vtkSmartPointer< T > New()
Create an instance of a VTK object.
std::vector< vtkHullFace > vtkHullFacesArray
vtkHullFace * GetFace(int faceId)
void ProducePolyData(vtkPolyData *pd)
Produce a vtkPolyData from the current polyhedron.
vtkHullFacesArray Faces
void AddFacePoint(vtkHullFace *face, vtkIdType ptId)
Add a point id defining the current face.
vtkHullPointsArray Points
vtkSmartPointer< vtkDoubleArray > Petals
ClipIntersectionStatus IntersectWithPlane(double origin[3], double normal[3], vtkIdType neiPtId)
The core geometric intersection operation.
FaceProcessingArray InProcessFaces
vtkVoronoiHull()
Constructor.
vtkIdType AddNewPoint(double x, double y, double z)
Methods to add and delete polyhedron points.
int EvaluateFace(vtkHullFace *face, int &startIdx, int &numKeptPts)
After a face clipping operation, characterize the face, and provide information for subsequent proces...
void UpdatePetals(double cf2)
void ComputeCircumFlower()
void SetPointFaces(vtkIdType pId, vtkIdType f0, vtkIdType f1, vtkIdType f2)
int IntersectFaceEdge(int faceId, int p0, int p1)
Given two point ids that form the edge of a polyhedron face, intersect the edge to produce a new inte...
DeletionStack DeletedFaces
std::vector< vtkHullPoint * > SortP
ClipIntersectionStatus Clip(vtkIdType neiPtId, const double neiPt[3])
Method to clip the current convex polyhedron/hull with a plane defined by a neighboring point.
int AddNewFace(vtkIdType npts, vtkIdType neiPtId)
Methods to create, modify, and delete polyhedron faces.
vtkVoronoiRandom01Range Bumper
std::vector< int > FaceScratchArray
void RebuildFace(int faceId, int startIdx, int numKeptPts)
Rebuild a convex, intersected face after a clipping operation.
FacePointsArray FacePoints
Internal data members.
bool InFlower(const double x[3])
void RebuildFacePoints(vtkHullFace *face, FaceScratchArray &idsBuffer)
After clipping, rebuild the face.
void Clear()
Empty out the polyhedron: clear memory but leave allocation intact.
int GetFacePoint(vtkHullFace *face, int ptNum)
Return the nth point id defining the current face.
void ProduceFacePolyData(vtkPolyData *pd, vtkHullFace *face)
Produce a vtkPolyData from the current polyhedron and one specified face.
void BumpNormal(int bumpNum, double normal[3], double bumpNormal[3])
void DeletePoint(vtkIdType ptId)
std::vector< int > FacePointsArray
vtkIdType AddNewPoint(double x[3])
DeletionStack DeletedPoints
void Initialize(vtkIdType genPtId, const double genPt[3], double bds[6])
Method to initiate the construction of the polyhedron.
void BumpOrigin(int bumpNum, double origin[3], double bumpOrigin[3])
FaceScratchArray FaceIdsBuffer
PointProcessingArray InProcessPoints
void MapPoints()
Used to produce debugging output (e.g., generate vtkPolyData).
vtkDoubleArray * GetPetals()
void DeleteFace(int faceId, int startIdx, int numKeptPts)
Delete a face from the polyhedron.
bool InCircumFlower(double r2)
std::vector< vtkHullEdgeTuple > InsertedEdgePointsArray
void AddNthFacePoint(vtkHullFace *face, int idx, vtkIdType ptId)
Add the nth point id defining the current face.
double GetCircumFlower2()
Methods to determine whether a point x[3] is within the Voronoi flower, or Voronoi circumflower.
InsertedEdgePointsArray InsertedEdgePoints
std::vector< vtkHullPoint > vtkHullPointsArray
void AllocatePointIds(int npts, vtkHullFace &face)
Internal memory operation to allocate space when adding new points (due to a reabuild) which define a...
bool operator!=(const vtkHullEdgeTuple &et) const
bool IsEdge(vtkIdType v0, vtkIdType v1) const
vtkIdType Faces[2]
vtkHullEdgeTuple(vtkIdType v0, vtkIdType v1)
bool operator==(const vtkHullEdgeTuple &et) const
Represent a face composing the polyhedron.
vtkIdType NeiId
void Replace(vtkIdType neiPtId)
vtkHullFace(vtkIdType neiId)
ProcessingStatus Status
Represent a point/vertex of the Voronoi polyhedron.
vtkHullPoint(double x, double y, double z)
ProcessingStatus Status
void Replace(double x, double y, double z)
vtkHullPoint(double x[3])
A homebrew stack with a preferred API.
Keep track of points and faces currently being operated on.
bool operator==(int faceId) const
std::function< void(vtkVoronoiHull &hull, int, int, int)> Function
Unique list of faces (by id) that require processing.
bool AddFace(vtkVoronoiHull *hull, int faceId)
Add points to be processed.
bool AddPoint(vtkVoronoiHull *hull, vtkHullPoint &point, int ptId)
@ Valid
Cell is in a good state.
int vtkIdType
Definition vtkType.h:363
#define VTK_FLOAT_MIN
Definition vtkType.h:199
#define VTK_FLOAT_MAX
Definition vtkType.h:200
ClipIntersectionStatus
Classes, structs, and typedefs in support of Voronoi processing.
ProcessingStatus
Clipping plane operations result in the dynamic deletion, modification, and addition of points and fa...