VTK  9.6.20260309
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
271 : PtId(-1)
272 , X{ 0, 0, 0 }
273 , NumClips(0)
274 , PruneTolerance(1.0e-13)
276 , RecomputePetals(true)
277 , CircumFlower2(0.0)
278 , MinRadius2(0.0)
279 , MaxRadius2(0.0)
280 {
281 // Preallocate some space
282 this->Points.reserve(256);
283 this->Faces.reserve(256);
284 this->FacePoints.reserve(2048);
285 this->InProcessPoints.reserve(256);
286 this->InProcessFaces.reserve(256);
287 this->DeletedPoints.Reserve(256);
288 this->DeletedFaces.Reserve(256);
289
290 // Each tile owns its own vtkDoubleArray
291 this->SortP.reserve(256);
293 this->Petals->SetNumberOfComponents(4); // x-y-z-R2
294 this->Petals->Allocate(256); // initial allocation
295 }
296
297 vtkVoronoiHull& operator=(const vtkVoronoiHull&) { return *this; }
298
304 void Initialize(vtkIdType genPtId, const double genPt[3], double bds[6]);
305
313 ClipIntersectionStatus Clip(vtkIdType neiPtId, const double neiPt[3]);
314
322 double GetCircumFlower2() { return this->CircumFlower2; }
323 bool InCircumFlower(double r2) // radius**2 of point from generator
324 {
325 // Only recompute the circumflower if necessary; that is, when
326 // a maximal point is eliminated by a polyhedral plane clip.
327 if (this->RecomputeCircumFlower)
328 {
329 this->ComputeCircumFlower();
330 }
332 }
333 bool InFlower(const double x[3]);
334 void UpdatePetals(double cf2);
336 {
337 if (this->RecomputePetals)
338 {
339 this->UpdatePetals(this->CircumFlower2);
340 }
341 return (this->Petals->GetNumberOfTuples() > 0 ? this->Petals : nullptr);
342 }
343
348 void MapPoints();
349
355
361
362 // Data members purposely left public for using classes to extract
363 // information.
364
365 // Information used to define the polyhedron- its generating point id and
366 // position, plus region classification. Indicate whether degenerate faces
367 // (i.e., those having ~zero area) can be deleted (i.e., pruned).
368 vtkIdType PtId; // Generating point id
369 double X[3]; // Generating point position
370 vtkIdType NumClips; // The total number of clip operations since Initialize()
371 double PruneTolerance; // Specify the prune tolerance
372
373 // Support normal joggle in the case of degeneracies.
375 void BumpNormal(int bumpNum, double normal[3], double bumpNormal[3]);
376 void BumpOrigin(int bumpNum, double origin[3], double bumpOrigin[3]);
377
378 // These data members represent the constructed polyhedron.
379 vtkIdType NumPts; // The number of valid points in the points array
380 vtkHullPointsArray Points; // Array of points defining this polyhedron
381 vtkIdType NumFaces; // The number of valid faces in the faces array
382 vtkHullFacesArray Faces; // A list of faces forming this polyhedron
383
390 {
391 std::vector<int> Stack;
392 DeletionStack() { this->Stack.reserve(256); }
393 void Push(int id) { this->Stack.emplace_back(id); }
394 int Pop()
395 {
396 vtkIdType id = this->Stack.back();
397 this->Stack.pop_back();
398 return id;
399 }
400 bool IsEmpty() { return (this->Stack.empty()); }
401 void Clear() { this->Stack.clear(); }
402 void Reserve(int sze) { this->Stack.reserve(sze); }
403 }; // DeletionStack
404
408
413 {
414 std::function<void(vtkVoronoiHull& hull, int, int, int)> Function;
418
419 // Constructor
420 FaceOperation(int faceId)
421 : FaceId(faceId)
422 {
423 }
424
425 // Support std::find()
426 bool operator==(int faceId) const { return (this->FaceId == faceId); }
427 };
428
435 struct FaceProcessingArray : public std::vector<FaceOperation>
436 {
437 bool AddFace(vtkVoronoiHull* hull, int faceId)
438 {
439 if (std::find(this->begin(), this->end(), faceId) == this->end())
440 {
441 FaceOperation& faceOp = this->emplace_back(faceId);
442 int numInts = hull->EvaluateFace(hull->GetFace(faceId), faceOp.StartIdx, faceOp.NumKeptPts);
443 if (numInts <= 0)
444 {
446 }
447 else if (numInts == 2)
448 {
450 }
451 else
452 {
453 return false;
454 }
455 }
456 return true; // successfully added the face for processing
457 }
458 };
459
463 struct PointProcessingArray : public std::vector<int>
464 {
465 // Add a non-degenerate point, and connected faces, for
466 // processing.
467 bool AddPoint(vtkVoronoiHull* hull, vtkHullPoint& point, int ptId)
468 {
469 this->emplace_back(ptId);
470 if (hull->InProcessFaces.AddFace(hull, point.Faces[0]) &&
471 hull->InProcessFaces.AddFace(hull, point.Faces[1]) &&
472 hull->InProcessFaces.AddFace(hull, point.Faces[2]))
473 {
474 return true;
475 }
476 else
477 {
478 return false;
479 }
480 }
481 };
482
486 vtkIdType AddNewPoint(double x, double y, double z);
487 vtkIdType AddNewPoint(double x[3]) { return this->AddNewPoint(x[0], x[1], x[2]); }
488 void DeletePoint(vtkIdType ptId);
489 void SetPointFaces(vtkIdType pId, vtkIdType f0, vtkIdType f1, vtkIdType f2);
490
494 int AddNewFace(vtkIdType npts, vtkIdType neiPtId);
495 vtkHullFace* GetFace(int faceId) { return &(this->Faces[faceId]); }
501 {
502 this->FacePoints[face->Points + face->InsertPos++] = ptId;
503 }
504
508 void AddNthFacePoint(vtkHullFace* face, int idx, vtkIdType ptId)
509 {
510 this->FacePoints[face->Points + idx] = ptId;
511 }
512
515 int GetFacePoint(vtkHullFace* face, int ptNum) { return this->FacePoints[face->Points + ptNum]; }
516
520 void RebuildFacePoints(vtkHullFace* face, FaceScratchArray& idsBuffer);
521
536 int EvaluateFace(vtkHullFace* face, int& startIdx, int& numKeptPts);
537
545 void DeleteFace(int faceId, int vtkNotUsed(startIdx), int vtkNotUsed(numKeptPts));
546
554 void RebuildFace(int faceId, int startIdx, int numKeptPts);
555
561 int IntersectFaceEdge(int faceId, int p0, int p1);
562
567 void AllocatePointIds(int npts, vtkHullFace& face);
568
569protected:
573 FacePointsArray FacePoints; // Point ids used to define faces
574 int MIN_POINTIDS_ALLOC = 10; // Minimum buffer allocation for face point ids
575
576 // Keep track of deleted points and faces so their memory can be reused. This acts as a
577 // sort of poor person's memory pool.
580
581 // Used to process and track faces and points affected by a plane clip operation.
582 PointProcessingArray InProcessPoints; // Points affected by current clip operation
583 FaceProcessingArray InProcessFaces; // Faces affected by current clip operation
584 InsertedEdgePointsArray InsertedEdgePoints; // New points generated on intersected edges
585 FaceScratchArray FaceIdsBuffer; // An internal buffer used to rebuild faces after clipping
586
587 // Indicate whether the Voronoi circumflower needs recomputing, and
588 // keep track of the current circumflower and related information.
589 void ComputeCircumFlower();
595 std::vector<vtkHullPoint*> SortP; // Points sorted on radius**2
596 vtkSmartPointer<vtkDoubleArray> Petals; // Flower petals w/ radii > shell radius
597
598 // Internal methods for processing clip operations of various types.
603 void Clear()
604 {
605 this->NumPts = 0;
606 this->Points.clear();
607 this->NumFaces = 0;
608 this->Faces.clear();
609 this->FacePoints.clear();
610 this->DeletedPoints.Clear();
611 this->DeletedFaces.Clear();
612 }
613
618 ClipIntersectionStatus IntersectWithPlane(double origin[3], double normal[3], vtkIdType neiPtId);
619}; // vtkVoronoiHull
620
621// In the following, inlined methods for performance
622
623//------------------------------------------------------------------------------
625{
626 this->Points[pId].Faces[0] = f0;
627 this->Points[pId].Faces[1] = f1;
628 this->Points[pId].Faces[2] = f2;
629}
630
631//------------------------------------------------------------------------------
632inline vtkIdType vtkVoronoiHull::AddNewPoint(double x, double y, double z)
633{
634 vtkIdType id;
635 // If there are no empty slots in the points array, then allocate a new point.
636 if (this->DeletedPoints.IsEmpty())
637 {
638 id = this->Points.size();
639 this->Points.emplace_back(x, y, z);
640 this->Points.back().R2 = vtkMath::Distance2BetweenPoints(this->Points.back().X, this->X);
641 }
642 else // otherwise, replace a previously deleted point and reuse its memory.
643 {
644 id = this->DeletedPoints.Pop();
645 this->Points[id].Replace(x, y, z);
646 this->Points[id].R2 = vtkMath::Distance2BetweenPoints(this->Points[id].X, this->X);
647 }
648
649 this->NumPts++;
650 return id;
651}
652
653//------------------------------------------------------------------------------
655{
656 // The circumflower may need recomputing if this point is an extreme point.
657 if ((4 * this->Points[ptId].R2) >= this->CircumFlower2)
658 {
659 this->RecomputeCircumFlower = true;
660 }
661
662 this->DeletedPoints.Push(ptId);
663 this->Points[ptId].Status = ProcessingStatus::Deleted;
664 this->NumPts--;
665}
666
667//------------------------------------------------------------------------------
669{
670 vtkIdType id;
671 // If there are no empty slots in the faces array, than allocate a new face.
672 if (this->DeletedFaces.IsEmpty())
673 {
674 id = this->Faces.size();
675 this->Faces.emplace_back(neiPtId);
676 this->AllocatePointIds(npts, this->Faces.back());
677 }
678 else
679 {
680 id = this->DeletedFaces.Pop();
681 this->Faces[id].Replace(neiPtId);
682 this->AllocatePointIds(npts, this->Faces[id]);
683 }
684
685 this->NumFaces++;
686 return id;
687}
688
689//------------------------------------------------------------------------------
690inline int vtkVoronoiHull::EvaluateFace(vtkHullFace* face, int& startIdx, int& numKeptPts)
691{
692 startIdx = numKeptPts = 0;
693 int npts = face->NumPts, numEdgeInts = 0;
694 int ip, p0, p1;
695 double val0, val1;
696
697 for (int i = 0; i < npts; ++i)
698 {
699 p0 = this->GetFacePoint(face, i);
700 ip = ((i + 1) == npts ? 0 : (i + 1));
701 p1 = this->GetFacePoint(face, ip);
702
703 val0 = this->Points[p0].Val;
704 val1 = this->Points[p1].Val;
705
706 if (val0 <= 0.0)
707 {
708 numKeptPts++;
709 }
710
711 if (val0 > 0 && val1 <= 0)
712 {
713 numEdgeInts++;
714 startIdx = i;
715 }
716 else if (val0 <= 0 && val1 > 0)
717 {
718 numEdgeInts++;
719 }
720 }
721
722 return numEdgeInts;
723}
724
725//------------------------------------------------------------------------------
727{
728 face.NumPts = npts;
729 face.InsertPos = 0;
730
731 // See if allocation is necessary. Otherwise use previous.
732 if (npts > face.AllocSize)
733 {
734 int size = (npts > MIN_POINTIDS_ALLOC ? npts : MIN_POINTIDS_ALLOC);
735 vtkIdType offset = this->FacePoints.size();
736 this->FacePoints.insert(this->FacePoints.end(), size, (-1));
737 face.AllocSize = size;
738 face.Points = offset;
739 }
740}
741
742//------------------------------------------------------------------------------
744{
745 // Make sure space has been allocated for rebuilt point ids.
746 int npts = static_cast<int>(idsBuffer.size());
747 if (npts > face->AllocSize)
748 {
749 // Need to realloc space for point ids.
750 this->AllocatePointIds(npts, *face);
751 }
752
753 // Copy buffer into face point ids.
754 for (face->InsertPos = 0; face->InsertPos < npts; ++face->InsertPos)
755 {
756 this->FacePoints[face->Points + face->InsertPos] = idsBuffer[face->InsertPos];
757 }
758 face->NumPts = npts;
759}
760
761//------------------------------------------------------------------------------
762inline void vtkVoronoiHull::DeleteFace(int faceId, int, int)
763{
764 this->DeletedFaces.Push(faceId);
765 this->Faces[faceId].Status = ProcessingStatus::Deleted;
766 this->NumFaces--;
767}
768
769//------------------------------------------------------------------------------
771{
772 // Compute the circumflower, and compute some info about
773 // the flower radii.
776
777 // Determine the circumflower and minimal sphere radius by
778 // checking against each of the flower petals.
779 for (const auto& pt : this->Points)
780 {
781 if (pt.Status == ProcessingStatus::Valid)
782 {
783 this->MinRadius2 = std::min(pt.R2, this->MinRadius2);
784 this->MaxRadius2 = std::max(pt.R2, this->MaxRadius2);
785 }
786 }
787 this->CircumFlower2 = (4.0 * this->MaxRadius2); // (2*(max petal radius))**2
788 this->RecomputeCircumFlower = false; // circumflower is up to date
789}
790
791//------------------------------------------------------------------------------
792inline bool vtkVoronoiHull::InFlower(const double x[3])
793{
794 // Check against the flower petals
795 for (const auto& pt : this->Points)
796 {
797 if (pt.Status == ProcessingStatus::Valid)
798 {
799 double r2 = vtkMath::Distance2BetweenPoints(x, pt.X);
800 if (r2 <= pt.R2)
801 {
802 return true;
803 }
804 }
805 } // for all valid points
806 // Point not in the flower since it's not in any petals
807 return false;
808}
809
810//------------------------------------------------------------------------------
812{
813 // Renumber the output points. Note that associated Faces should use the
814 // PtMap id to ensure the the point connectivity ids are contiguous.
815 vtkIdType id = 0;
816 for (auto& pitr : this->Points)
817 {
818 if (pitr.Status == ProcessingStatus::Valid)
819 {
820 pitr.PtMap = id++;
821 }
822 }
823}
824
825//------------------------------------------------------------------------------
826// Update the flower petals which are passed off to the locator.
827// Only petals which exend past the minimal radius of the shell
828// request are added to the list of petals. It is presumed that
829// UpdateCircumFlower() has been invoked previously.
830inline void vtkVoronoiHull::UpdatePetals(double cf2)
831{
832 // If the radii of the flower spheres (petals) is highly variable (which
833 // occurs when the spacing of points is highly variable), then there is
834 // likely a lot of empty search space. Only add flower petals which extend
835 // past the outer shell request boundary. These petals are used to further
836 // limit the point search space.
837 this->Petals->Reset();
838 this->RecomputePetals = false; // petals will be updated in the following
839
840 constexpr double SphereRatio = 2.0, SphereRatio2 = SphereRatio * SphereRatio;
841 if (this->MinRadius2 > 0 && ((this->MaxRadius2 / this->MinRadius2) < SphereRatio2))
842 {
843 return; // it's not worth using the petals
844 }
845
846 // Empirically determined
847 constexpr double LargeSphereRatio = 0.25;
848 int maxLargeSpheres = LargeSphereRatio * this->NumPts;
849
850 this->SortP.clear();
851 double minR2 = VTK_FLOAT_MAX, maxR2 = VTK_FLOAT_MIN;
852 for (auto& pt : this->Points)
853 {
854 if (pt.Status == ProcessingStatus::Valid)
855 {
856 // (2*R)**2 >= shell request radius**2
857 if ((4 * pt.R2) >= cf2)
858 {
859 minR2 = std::min(pt.R2, minR2);
860 maxR2 = std::max(pt.R2, maxR2);
861 this->SortP.emplace_back(&pt);
862 }
863 }
864 }
865
866 if (static_cast<int>(this->SortP.size()) > maxLargeSpheres || (maxR2 / minR2) < SphereRatio2)
867 {
868 return; // it's not worth using the petals
869 }
870 else // sort from large spheres to small
871 {
872 std::sort(this->SortP.begin(), this->SortP.end(),
873 [](const vtkHullPoint* p0, const vtkHullPoint* p1) { return (p0->R2 > p1->R2); });
874 for (auto& pt : this->SortP)
875 {
876 this->Petals->InsertNextTuple4(pt->X[0], pt->X[1], pt->X[2], pt->R2);
877 }
878 } // it's worth using large sphere culling
879} // UpdatePetals
880
881VTK_ABI_NAMESPACE_END
882#endif
883// 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
vtkVoronoiHull(const vtkVoronoiHull &)
vtkSMPTools require a proper operator= and copy constructor.
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...
vtkVoronoiHull & operator=(const vtkVoronoiHull &)
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...