VTK  9.6.20260212
vtkVoronoiCore.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
18
19#ifndef vtkVoronoiCore_h
20#define vtkVoronoiCore_h
21
22#include "vtkAlgorithm.h" // check abort status if embedded in filter
23#include "vtkIntArray.h" // for representing segmented region ids
24#include "vtkSMPTools.h" // SMP parallel processing
25
26#include <algorithm> // for std::sort
27#include <iostream>
28#include <random> // random generation of colors
29#include <vector>
30
31VTK_ABI_NAMESPACE_BEGIN
32
45
52{
53 BACKWARD_SPOKE = 0, // Bit 0: Backward spoke
54 FORWARD_SPOKE = 1, // Bit 0: Forward spoke
55 REGION_BOUNDARY = 2, // Bit 1: Region boundary spoke
56 DOMAIN_BOUNDARY = 4, // Bit 2: Domain boundary spoke
57 PRUNED = 8, // Bit 3: Spoke is pruned (deleted)
58};
59
64{
65 vtkIdType NeiId; // Id of the wheel that the spoke is connected to (wheelId,NeiId)
66 unsigned char Classification; // Indicate the classification of this spoke
68 : NeiId(-1)
70 {
71 }
72 vtkVoronoiSpoke(vtkIdType neiId, unsigned char classification)
73 : NeiId(neiId)
74 , Classification(classification)
75 {
76 }
77}; // vtkVoronoiSpoke
78
86using vtkVoronoiSpokesType = std::vector<vtkVoronoiSpoke>;
87using vtkVoronoiSpokesIterator = vtkVoronoiSpokesType::iterator;
88using vtkVoronoiWheelsType = std::vector<vtkIdType>;
89
90// Gather spokes into a wheel. Define some basic operators. Note that every
91// wheel is associated with an input (tile/hull generating) point. So access
92// to the wheel and its associated spokes is via point id.
94{
95 vtkVoronoiWheelsType& Wheels; // The composited array of wheels
96 vtkVoronoiSpokesType& Spokes; // The composited array of spokes
97 vtkIdType Id; // The associated point/tile id: with wheelId == pointId
98 int NumSpokes; // The number of emanating spokes
99 vtkVoronoiSpoke* LocalSpokes; // The array of spokes connected to this wheel
100
101 // Default instantiation.
103 : Wheels(wheels)
104 , Spokes(spokes)
105 , Id(-1)
106 , NumSpokes(0)
107 , LocalSpokes(nullptr)
108 {
109 }
110
111 // Setup the wheel for queries: an efficient form that does not require
112 // repeated wheel instantiation.
114 {
115 this->Id = id;
116 this->NumSpokes = numSpokes = (this->Wheels[id + 1] - this->Wheels[id]);
117 this->LocalSpokes = &(this->Spokes[this->Wheels[id]]);
118 return this->LocalSpokes;
119 }
120}; // vtkVoronoiWheel
121
130{
131 vtkVoronoiWheelsType Wheels; // Wheel/spokes data structure: offset array into spokes
132 vtkVoronoiSpokesType Spokes; // Spokes / edges with classification
133
134 void Initialize(vtkIdType numWheels, vtkIdType numSpokes);
135 vtkVoronoiWheelsType& GetWheels() { return this->Wheels; }
136 vtkVoronoiSpokesType& GetSpokes() { return this->Spokes; }
137 vtkIdType GetNumberOfWheels() { return (this->Wheels.size() - 1); }
138 vtkIdType GetNumberOfSpokes() { return this->Spokes.size(); }
140 vtkIdType GetNumberOfSpokes(vtkIdType ptId); // #spokes for a specified wheel
142 vtkIdType GetWheelOffset(vtkIdType ptId) { return this->Wheels[ptId]; }
143 static void CountFaces(const vtkVoronoiSpoke* spokes, int numSpokes, int& numDomainBoundaryFaces,
144 int& numRegionBoundaryFaces, int& numForwardFaces);
145
153 bool Validate();
154
159 {
162
164 : Graph(graph)
165 , NumInvalid(0)
166 {
167 }
168
169 // Keep track whether threads are non-degenerate.
171
172 // vtkSMPTools threaded interface
174 void operator()(vtkIdType wheelId, vtkIdType endWheelId);
175 void Reduce();
176 }; // ValidateAdjacencyGraph
177}; // vtkVoronoiAdjacencyGraph
178
183struct vtkVoronoiHullVertex // 3D hull vertices
184{
185 double X[3];
186 vtkVoronoiHullVertex(double x, double y, double z)
187 : X{ x, y, z }
188 {
189 }
190 vtkVoronoiHullVertex(const double x[3])
191 : X{ x[0], x[1], x[2] }
192 {
193 }
194};
195using vtkVoronoiHullVertexType = std::vector<vtkVoronoiHullVertex>;
196
197struct vtkVoronoiTileVertex // 2D tile vertices
198{
199 double X[2];
200 vtkVoronoiTileVertex(double x, double y)
201 : X{ x, y }
202 {
203 }
204 vtkVoronoiTileVertex(const double x[2])
205 : X{ x[0], x[1] }
206 {
207 }
208};
209using vtkVoronoiTileVertexType = std::vector<vtkVoronoiTileVertex>;
210
222{
227 std::array<vtkIdType, 4> Ids;
228
233 : Ids{ 0 }
234 {
235 }
236
241 : Ids{ p0, p1, p2, ptId }
242 {
243 std::sort(this->Ids.data(), this->Ids.data() + this->Ids.size());
244 }
245
248 vtkVoronoiTopoCoord3D(const vtkVoronoiTopoCoord3D& tt) { this->Ids = tt.Ids; }
249
255 bool operator<(const vtkVoronoiTopoCoord3D& tuple) const { return this->Ids < tuple.Ids; }
256};
257using vtkVoronoiTopoCoords3DType = std::vector<vtkVoronoiTopoCoord3D>;
258
260{
265 std::array<vtkIdType, 3> Ids;
266
271 : Ids{ 0 }
272 {
273 }
274
279 : Ids{ p0, p1, ptId }
280 {
281 std::sort(this->Ids.data(), this->Ids.data() + this->Ids.size());
282 }
283
286 vtkVoronoiTopoCoord2D(const vtkVoronoiTopoCoord2D& tt) { this->Ids = tt.Ids; }
287
293 bool operator<(const vtkVoronoiTopoCoord2D& tuple) const { return this->Ids < tuple.Ids; }
294};
295using vtkVoronoiTopoCoords2DType = std::vector<vtkVoronoiTopoCoord2D>;
296
305{
306 vtkIdType PtId; // the id of the hull vertex
307
309 : PtId(-1)
310 {
311 }
312 bool operator!=(const vtkVoronoiMergeTuple3D& mt) const { return (this->Ids != mt.Ids); }
313}; // vtkVoronoiMergeTuple3D
314
316{
317 vtkIdType PtId; // the id of the tile vertex
318
320 : PtId(-1)
321 {
322 }
323 bool operator!=(const vtkVoronoiMergeTuple2D& mt) const { return (this->Ids != mt.Ids); }
324}; // vtkVoronoiMergeTuple2D
325
332using vtkMergeTupleOffsets = std::vector<vtkIdType>; // offsets into merged tuples
333using vtkMergeTuples3DType = std::vector<vtkVoronoiMergeTuple3D>;
334using vtkMergeTuples2DType = std::vector<vtkVoronoiMergeTuple2D>;
335
341using vtkMergeMapType = std::vector<vtkIdType>;
342
346using vtkVoronoiCellConnType = std::vector<vtkIdType>;
347
353{
354 vtkIdType Num; // Number of total items (e.g., points) to process
355 vtkIdType BatchSize; // The desired batch size (clamped by Num)
356 vtkIdType NumBatches; // The total number of batches to process
358 : Num(num)
359 , BatchSize(batchSize)
360 {
361 this->NumBatches = static_cast<vtkIdType>(ceil(static_cast<double>(num) / batchSize));
362 }
363 vtkIdType GetNumberOfBatches() const { return this->NumBatches; }
364 vtkIdType GetBatchItemRange(vtkIdType batchNum, vtkIdType& startId, vtkIdType& endId) const
365 {
366 startId = batchNum * this->BatchSize;
367 endId = startId + this->BatchSize;
368 endId = (endId > this->Num ? this->Num : endId);
369 return (endId - startId);
370 }
371}; // vtkVoronoiBatchManager
372
377using vtkBatchIdsType = std::vector<vtkIdType>;
378
379// Convenience function: convert input labels/region ids/scalars to signed int.
380// The Voronoi classes expect signed int region labels.
382{
384 rIds->SetNumberOfTuples(inScalars->GetNumberOfTuples());
385 rIds->DeepCopy(inScalars);
386
387 return rIds;
388}
389
395{
399
401 : Filter(filter)
402 {
403 this->IsFirst = vtkSMPTools::GetSingleThread();
404 this->CheckAbortInterval = std::min((end - start) / 10 + 1, (vtkIdType)1000);
405 }
406
408 {
409 if (this->Filter && this->IsFirst && !(id % this->CheckAbortInterval))
410 {
411 this->Filter->CheckAbort();
412 return (this->Filter->GetAbortOutput() ? true : false);
413 }
414 return false;
415 }
416}; // vtkVoronoiAbortCheck
417
418// Use system <random> - create a simple convenience class. This generates
419// random color indices [0,64).
421{
422 std::mt19937 RNG;
423 std::uniform_int_distribution<int> Dist;
424 vtkVoronoiRandomColors() { this->Dist.param(typename decltype(this->Dist)::param_type(0, 64)); }
425 void Seed(vtkIdType s) { this->RNG.seed(s); }
426 vtkIdType Next() { return this->Dist(RNG); }
427};
428
429// Use system <random> - create a simple convenience class. This generates
430// random real values [0,1).
432{
433 std::mt19937 RNG;
434 std::uniform_real_distribution<double> Dist;
436 {
437 this->Dist.param(typename decltype(this->Dist)::param_type(0.0, 1.0));
438 }
439 void Seed(vtkIdType s) { this->RNG.seed(s); }
440 double Next() { return this->Dist(RNG); }
441};
442
443// A convenience class and methods to randomly perturb (joggle or jitter)
444// point positions. Such jittering (even if very small) significantly
445// improves the numerical stability of Voronoi and Delaunay computations.
446// Implementation note: these methods might be better added to vtkMath since
447// they can be used by other classes. Once they are demonstrated to be stable,
448// they may be moved.
450{
451 // Joggle a single point at input position xIn to produce the output position
452 // xOut (xIn and xOut may be computed in-place). The radius is the allowable
453 // range of joggle in the sphere. A sequence is provided, assumed properly
454 // initialized, to produce random (0,1) values. Note that if this method is
455 // invoked in a thread, separate sequence instantiations (one per thread)
456 // should be provided.
457 static void JoggleXYZ(
458 double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range& sequence)
459 {
460 double cosphi = 1 - 2 * sequence.Next();
461 double sinphi = sqrt(1 - cosphi * cosphi);
462 double rho = radius * pow(sequence.Next(), 0.33333333);
463 double R = rho * sinphi;
464 double theta = 2.0 * vtkMath::Pi() * sequence.Next();
465 xOut[0] = xIn[0] + R * cos(theta);
466 xOut[1] = xIn[1] + R * sin(theta);
467 xOut[2] = xIn[2] + rho * cosphi;
468 }
469
470 // Joggle a single point at input position xIn to produce the output
471 // position xOut (xIn and xOut may be computed in-place). The radius is the
472 // allowable range of joggle in the circle in the x-y plane. A sequence is
473 // provided, assumed properly initialized, to produce random (0,1) values.
474 // Note that if this method is invoked in a thread, separate sequence
475 // instantiations (one per thread) should be provided.
476 static void JoggleXY(
477 double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range& sequence)
478 {
479 double R = radius * sequence.Next();
480 double theta = 2.0 * vtkMath::Pi() * sequence.Next();
481 xOut[0] = xIn[0] + R * cos(theta);
482 xOut[1] = xIn[1] + R * sin(theta);
483 xOut[2] = xIn[2];
484 }
485
486 // Joggle a single point at input position xIn to produce the output
487 // position xOut (xIn and xOut may be computed in-place). The radius is the
488 // allowable range of joggle in the circle in the x-z plane. A sequence is
489 // provided, assumed properly initialized, to produce random (0,1) values.
490 // Note that if this method is invoked in a thread, separate sequence
491 // instantiations (one per thread) should be provided.
492 static void JoggleXZ(
493 double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range& sequence)
494 {
495 double R = radius * sequence.Next();
496 double theta = 2.0 * vtkMath::Pi() * sequence.Next();
497 xOut[0] = xIn[0] + R * cos(theta);
498 xOut[1] = xIn[1];
499 xOut[2] = xIn[2] + R * sin(theta);
500 }
501
502 // Joggle a single point at input position xIn to produce the output
503 // position xOut (xIn and xOut may be computed in-place). The radius is the
504 // allowable range of joggle in the circle in the y-z plane. A sequence is
505 // provided, assumed properly initialized, to produce random (0,1) values.
506 // Note that if this method is invoked in a thread, separate sequence
507 // instantiations (one per thread) should be provided.
508 static void JoggleYZ(
509 double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range& sequence)
510 {
511 double R = radius * sequence.Next();
512 double theta = 2.0 * vtkMath::Pi() * sequence.Next();
513 xOut[0] = xIn[0];
514 xOut[1] = xIn[1] + R * cos(theta);
515 xOut[2] = xIn[2] + R * sin(theta);
516 }
517
518 // Joggle a 3D vector. Specify an input normalized vector (vecIn), a angle
519 // of rotation theta (in radians, 0<theta<Pi/2), and a random
520 // sequence. Produce on output a normalized vector (vecOut) -- vecIn and
521 // vecOut may be computed in place. Note that if this method is invoked in
522 // a thread, separate sequence instantiations (one per thread) should be
523 // provided. The approach used here is to define a disk located at the tip
524 // of the (normalized) vecIn (the disc is normal to vecIn), and randomly
525 // select a vector from the edge/along the perimeter of the disc.
526 static void JoggleNormal(
527 double vecIn[3], double vecOut[3], double theta, vtkVoronoiRandom01Range& sequence)
528 {
529 // Find a vector orthogonal to the input vector.
530 double perp[3];
531 int iMax = (std::fabs(vecIn[0]) > std::fabs(vecIn[1]) ? 0 : 1);
532 iMax = (std::fabs(vecIn[iMax]) > std::fabs(vecIn[2]) ? iMax : 2);
533 perp[(iMax + 1) % 3] = 1.0;
534 perp[(iMax + 2) % 3] = 0.0;
535 perp[iMax] = -(vecIn[(iMax + 1) % 3] / vecIn[iMax]);
536 vtkMath::Normalize(perp);
537
538 // Now find third orthogonal vector. Since vecIn and perp are normalized and
539 // orthogonal, the resultant cross also a unit vector.
540 double cross[3];
541 vtkMath::Cross(vecIn, perp, cross);
542
543 // Compute a random vector
544 double radius = theta / (2.0 * vtkMath::Pi()); // small angle approximation
545 double t = 2.0 * vtkMath::Pi() * sequence.Next();
546 double x = radius * cos(t);
547 double y = radius * sin(t);
548
549 // Construct a vector along the disc/cone edge
550 vecOut[0] = vecIn[0] + x * perp[0] + y * cross[0];
551 vecOut[1] = vecIn[1] + x * perp[1] + y * cross[1];
552 vecOut[2] = vecIn[2] + x * perp[2] + y * cross[2];
553 vtkMath::Normalize(vecOut);
554 }
555}; // vtkVoronoiJoggle
556
557VTK_ABI_NAMESPACE_END
558#include "vtkVoronoiCore.txx"
559
560#endif
561// VTK-HeaderTest-Exclude: vtkVoronoiCore.h
RealT mt
Definition PyrC2Basis.h:39
vtkIdType GetNumberOfTuples() const
Get the number of complete tuples (a component group) in the array.
bool CheckAbort()
Checks to see if this filter should abort.
virtual bool GetAbortOutput()
Set/Get an internal variable used to communicate between the algorithm and executive.
static constexpr double Pi()
A mathematical constant.
Definition vtkMath.h:227
static float Normalize(float v[3])
Normalize (in place) a 3-vector.
Definition vtkMath.h:1972
static void Cross(VectorT1 &&a, VectorT2 &&b, VectorT3 &c)
Cross product of two 3-vectors.
Definition vtkMath.h:2079
Allocate and hold a VTK object.
Definition vtkNew.h:167
Thread local storage for VTK objects.
static bool GetSingleThread()
Returns true if the given thread is specified thread for single scope.
Hold a reference to a vtkObjectBase instance.
vtkVoronoiAbortCheck(vtkIdType start, vtkIdType end, vtkAlgorithm *filter)
bool operator()(vtkIdType id)
vtkAlgorithm * Filter
ValidateAdjacencyGraph(vtkVoronoiAdjacencyGraph &graph)
void operator()(vtkIdType wheelId, vtkIdType endWheelId)
vtkSMPThreadLocal< vtkIdType > ThreadNumInvalid
The adjacency graph, a collection of wheels and spokes, is a topological construct that connects Voro...
vtkVoronoiWheelsType Wheels
vtkIdType GetNumberOfSpokes(vtkIdType ptId)
vtkVoronoiSpoke * GetSpokes(vtkIdType ptId, vtkIdType &numSpokes)
bool Validate()
Return true if the graph meets the conditions necessary to form a valid Voronoi tessellation.
vtkVoronoiWheelsType & GetWheels()
static void CountFaces(const vtkVoronoiSpoke *spokes, int numSpokes, int &numDomainBoundaryFaces, int &numRegionBoundaryFaces, int &numForwardFaces)
void Initialize(vtkIdType numWheels, vtkIdType numSpokes)
vtkVoronoiSpokesType Spokes
bool IsSpoke(vtkIdType pt0, vtkIdType pt1)
vtkIdType GetWheelOffset(vtkIdType ptId)
vtkVoronoiSpokesType & GetSpokes()
vtkIdType GetBatchItemRange(vtkIdType batchNum, vtkIdType &startId, vtkIdType &endId) const
vtkVoronoiBatchManager(vtkIdType num, vtkIdType batchSize)
vtkIdType GetNumberOfBatches() const
vtkVoronoiHullVertex(const double x[3])
vtkVoronoiHullVertex(double x, double y, double z)
static void JoggleNormal(double vecIn[3], double vecOut[3], double theta, vtkVoronoiRandom01Range &sequence)
static void JoggleXY(double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range &sequence)
static void JoggleXZ(double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range &sequence)
static void JoggleYZ(double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range &sequence)
static void JoggleXYZ(double xIn[3], double xOut[3], double radius, vtkVoronoiRandom01Range &sequence)
bool operator!=(const vtkVoronoiMergeTuple2D &mt) const
bool operator!=(const vtkVoronoiMergeTuple3D &mt) const
std::uniform_real_distribution< double > Dist
void Seed(vtkIdType s)
std::uniform_int_distribution< int > Dist
void Seed(vtkIdType s)
Typedefs and classes in support of the adjacency graph.
unsigned char Classification
vtkVoronoiSpoke(vtkIdType neiId, unsigned char classification)
vtkVoronoiTileVertex(double x, double y)
vtkVoronoiTileVertex(const double x[2])
vtkVoronoiTopoCoord2D(vtkIdType p0, vtkIdType p1, vtkIdType ptId)
Define with the N+1 point generators: the N generators producing the hull vertex, plus the current po...
vtkVoronoiTopoCoord2D()
Various flavors of constructors.
vtkVoronoiTopoCoord2D(const vtkVoronoiTopoCoord2D &tt)
Copy constructor assumes that tuple ids are already sorted.
bool operator<(const vtkVoronoiTopoCoord2D &tuple) const
Operator< used to support a subsequent sort operation of the n-tuples (used for uniquely identifying ...
std::array< vtkIdType, 3 > Ids
Points defining a topological coord tuple / Delaunay simplex.
vtkVoronoiTopoCoord3D()
Various flavors of constructors.
bool operator<(const vtkVoronoiTopoCoord3D &tuple) const
Operator< used to support a subsequent sort operation of the n-tuples (used for uniquely identifying ...
std::array< vtkIdType, 4 > Ids
Points defining a topological coord tuple / Delaunay simplex.
vtkVoronoiTopoCoord3D(vtkIdType p0, vtkIdType p1, vtkIdType p2, vtkIdType ptId)
Define with the N+1 point generators: the N generators producing the hull vertex, plus the current po...
vtkVoronoiTopoCoord3D(const vtkVoronoiTopoCoord3D &tt)
Copy constructor assumes that tuple ids are already sorted.
vtkVoronoiSpoke * LocalSpokes
vtkVoronoiSpokesType & Spokes
vtkVoronoiWheel(vtkVoronoiWheelsType &wheels, vtkVoronoiSpokesType &spokes)
vtkVoronoiSpoke * Initialize(vtkIdType id, int &numSpokes)
vtkVoronoiWheelsType & Wheels
#define vtkDataArray
int vtkIdType
Definition vtkType.h:354
std::vector< vtkVoronoiTopoCoord2D > vtkVoronoiTopoCoords2DType
std::vector< vtkIdType > vtkVoronoiCellConnType
Convenience type for representing cell connectivity during compositing.
std::vector< vtkIdType > vtkBatchIdsType
Keep track of batches of generating points.
std::vector< vtkIdType > vtkMergeMapType
When merging points, the merge map is a vector that maps global tile/hull vertex ids (which contain d...
std::vector< vtkIdType > vtkVoronoiWheelsType
std::vector< vtkVoronoiSpoke > vtkVoronoiSpokesType
The vtkVoronoiWheelsType vector is used to keep track of the number of spokes (and equivalently,...
ClipIntersectionStatus
Classes, structs, and typedefs in support of Voronoi processing.
vtkSmartPointer< vtkIntArray > ConvertRegionLabels(vtkDataArray *inScalars)
vtkVoronoiSpokesType::iterator vtkVoronoiSpokesIterator
std::vector< vtkVoronoiHullVertex > vtkVoronoiHullVertexType
std::vector< vtkVoronoiMergeTuple3D > vtkMergeTuples3DType
vtkSpokeClassification
Classification for Voronoi spokes (and associated faces).
@ BACKWARD_SPOKE
@ PRUNED
@ DOMAIN_BOUNDARY
@ FORWARD_SPOKE
@ REGION_BOUNDARY
std::vector< vtkVoronoiMergeTuple2D > vtkMergeTuples2DType
std::vector< vtkVoronoiTopoCoord3D > vtkVoronoiTopoCoords3DType
std::vector< vtkIdType > vtkMergeTupleOffsets
Global tile/hull vertices, with duplicates, that are assigned a global id (if point merging is perfor...
std::vector< vtkVoronoiTileVertex > vtkVoronoiTileVertexType
#define vtkAlgorithm