TTK
Loading...
Searching...
No Matches
ttkTriangulationFactory.cpp
Go to the documentation of this file.
2
3#include <Triangulation.h>
4#include <ttkUtils.h>
5
6#include <vtkCallbackCommand.h>
7#include <vtkCellData.h>
8#include <vtkCellTypes.h>
9#include <vtkCommand.h>
10#include <vtkImageData.h>
11#include <vtkPointData.h>
12#include <vtkPolyData.h>
13#include <vtkUnstructuredGrid.h>
14#include <vtkVersionMacros.h>
15
16static vtkCellArray *GetCells(vtkDataSet *dataSet) {
17 switch(dataSet->GetDataObjectType()) {
18 case VTK_UNSTRUCTURED_GRID: {
19 auto dataSetAsUG = static_cast<vtkUnstructuredGrid *>(dataSet);
20 return dataSetAsUG->GetCells();
21 }
22 case VTK_POLY_DATA: {
23 auto dataSetAsPD = static_cast<vtkPolyData *>(dataSet);
24 return dataSetAsPD->GetNumberOfPolys() > 0 ? dataSetAsPD->GetPolys()
25 : dataSetAsPD->GetNumberOfLines() > 0 ? dataSetAsPD->GetLines()
26 : dataSetAsPD->GetVerts();
27 }
28 }
29 return nullptr;
30}
31
32static int checkCellTypes(vtkPointSet *object) {
33
34 size_t nTypes = 0;
35
36#if VTK_VERSION_NUMBER >= VTK_VERSION_CHECK(9, 2, 0)
37 if(object->GetDataObjectType() == VTK_UNSTRUCTURED_GRID) {
38 auto objectAsUG = vtkUnstructuredGrid::SafeDownCast(object);
39 auto distinctCellTypes = objectAsUG->GetDistinctCellTypesArray();
40 nTypes = distinctCellTypes->GetNumberOfTuples();
41 } else
42#endif
43 {
44 auto cellTypes = vtkSmartPointer<vtkCellTypes>::New();
45 object->GetCellTypes(cellTypes);
46 nTypes = cellTypes->GetNumberOfTypes();
47 }
48
49 // if cells are empty
50 if(nTypes == 0)
51 return 1; // no error
52
53 // if cells are not homogeneous
54 if(nTypes > 1)
55 return -1;
56
57 // if cells are not simplices
58 if(nTypes == 1) {
59 const auto &cellType = object->GetCellType(0);
60 if(cellType != VTK_VERTEX && cellType != VTK_LINE
61 && cellType != VTK_TRIANGLE && cellType != VTK_TETRA)
62 return -2;
63 }
64
65 return 1;
66}
67
68struct ttkOnDeleteCommand : public vtkCommand {
70 vtkObject *observee;
71
73 return new ttkOnDeleteCommand;
74 }
75 vtkTypeMacro(ttkOnDeleteCommand, vtkCommand);
76
77 void Init(vtkDataSet *dataSet) {
78 this->key = ttkTriangulationFactory::GetKey(dataSet);
79
80 if(dataSet->IsA("vtkPointSet"))
81 this->observee = static_cast<vtkObject *>(GetCells(dataSet));
82 else
83 this->observee = static_cast<vtkObject *>(dataSet);
84
85 this->observee->AddObserver(vtkCommand::DeleteEvent, this, 1);
86 }
87
88 void Execute(vtkObject *,
89 unsigned long ttkNotUsed(eventId),
90 void *ttkNotUsed(callData)) override {
91 if(this->observee)
92 this->observee->RemoveObserver(this);
93
94 auto instance = &ttkTriangulationFactory::Instance;
95
96 if(instance->registry.empty()) {
97 return;
98 }
99
100 auto it = instance->registry.find(this->key);
101 if(it != instance->registry.end()) {
102 instance->registry.erase(it);
103 instance->printMsg("Triangulation Deleted", ttk::debug::Priority::DETAIL);
104 instance->printMsg("# Registered Triangulations: "
105 + std::to_string(instance->registry.size()),
107 }
108 }
109};
110
112 ttk::Triangulation *triangulation_)
113 : triangulation(triangulation_), owner(dataSet) {
114 auto cells = GetCells(dataSet);
115 if(cells)
116 this->cellModTime = cells->GetMTime();
117
118 if(dataSet->IsA("vtkImageData")) {
119 auto image = static_cast<vtkImageData *>(dataSet);
120 image->GetExtent(this->extent);
121 image->GetOrigin(this->origin);
122 image->GetSpacing(this->spacing);
123 image->GetDimensions(this->dimensions);
124 }
125
127 onDelete->Init(dataSet);
128}
129
130bool RegistryValue::isValid(vtkDataSet *dataSet) const {
131 auto cells = GetCells(dataSet);
132 if(cells)
133 return this->cellModTime == cells->GetMTime();
134
135 if(dataSet->IsA("vtkImageData")) {
136 auto image = static_cast<vtkImageData *>(dataSet);
137
138 int extent_[6];
139 double origin_[3];
140 double spacing_[3];
141 int dimensions_[3];
142
143 image->GetExtent(extent_);
144 image->GetOrigin(origin_);
145 image->GetSpacing(spacing_);
146 image->GetDimensions(dimensions_);
147
148#ifdef TTK_ENABLE_MPI
149 bool isValid = triangulation->getIsMPIValid();
150#else
151 bool isValid = true;
152#endif
153 for(int i = 0; i < 6; i++)
154 if(this->extent[i] != extent_[i])
155 isValid = false;
156 for(int i = 0; i < 3; i++)
157 if(this->origin[i] != origin_[i] || this->spacing[i] != spacing_[i]
158 || this->dimensions[i] != dimensions_[i])
159 isValid = false;
160
161 return isValid;
162 }
163
164 return false;
165}
166
167ttkTriangulationFactory::ttkTriangulationFactory() {
168 this->setDebugMsgPrefix("TriangulationFactory");
169}
170
172 ttkTriangulationFactory::CreateImplicitTriangulation(vtkImageData *image) {
173 ttk::Timer timer;
174 this->printMsg("Initializing Implicit Triangulation", 0, 0,
176
177 auto triangulation = std::make_unique<ttk::Triangulation>();
178
179 int extent[6];
180 image->GetExtent(extent);
181
182 double origin[3];
183 image->GetOrigin(origin);
184
185 double spacing[3];
186 image->GetSpacing(spacing);
187
188 // 1D (not tested)
189 if(!spacing[1] && !spacing[2])
190 spacing[1] = 1;
191 if(!spacing[2]) // 2D (tested)
192 spacing[2] = 1;
193
194 int dimensions[3];
195 image->GetDimensions(dimensions);
196
197 double firstPoint[3];
198 firstPoint[0] = origin[0] + extent[0] * spacing[0];
199 firstPoint[1] = origin[1] + extent[2] * spacing[1];
200 firstPoint[2] = origin[2] + extent[4] * spacing[2];
201
202 triangulation->setInputGrid(firstPoint[0], firstPoint[1], firstPoint[2],
203 spacing[0], spacing[1], spacing[2], dimensions[0],
204 dimensions[1], dimensions[2]);
205
206 this->printMsg("Initializing Implicit Triangulation", 1,
209
210 return triangulation;
211}
212
214 ttkTriangulationFactory::CreateExplicitTriangulation(vtkPointSet *pointSet) {
215 ttk::Timer timer;
216
217 auto points = pointSet->GetPoints();
218 if(!points) {
219 this->printErr("DataSet has uninitialized `vtkPoints`.");
220 return nullptr;
221 }
222
223 auto cells = GetCells(pointSet);
224 if(!cells) {
225 this->printErr("DataSet has uninitialized `vtkCellArray`.");
226 return nullptr;
227 }
228
229 auto triangulation = std::make_unique<ttk::Triangulation>();
230 int const hasIndexArray
231 = pointSet->GetPointData()->HasArray(ttk::compactTriangulationIndex);
232
233 if(hasIndexArray) {
234 this->printMsg("Initializing Compact Triangulation", 0, 0,
236 } else {
237 this->printMsg("Initializing Explicit Triangulation", 0, 0,
239 }
240
241 // Points
242 {
243 auto pointDataType = points->GetDataType();
244 if(pointDataType != VTK_FLOAT && pointDataType != VTK_DOUBLE) {
245 this->printErr("Unable to initialize 'ttk::Triangulation' for point "
246 "precision other than 'float' or 'double'.");
247 return {};
248 }
249
250 void *pointDataArray = ttkUtils::GetVoidPointer(points);
251 if(hasIndexArray) {
252 vtkAbstractArray *indexArray = pointSet->GetPointData()->GetAbstractArray(
254 triangulation->setStellarInputPoints(
255 points->GetNumberOfPoints(), pointDataArray,
256 (int *)indexArray->GetVoidPointer(0), pointDataType == VTK_DOUBLE);
257 } else {
258 triangulation->setInputPoints(points->GetNumberOfPoints(), pointDataArray,
259 pointDataType == VTK_DOUBLE);
260 }
261 }
262
263 // check if cell types are simplices
264 int const cellTypeStatus = checkCellTypes(pointSet);
265 if(cellTypeStatus == -1) {
266 this->printWrn("Inhomogeneous cell dimensions detected.");
267 this->printWrn(
268 "Consider using `ttkExtract` to extract cells of a given dimension.");
269 return {};
270 } else if(cellTypeStatus == -2) {
271 this->printWrn("Cells are not simplices.");
272 this->printWrn("Consider using `vtkTetrahedralize` in pre-processing.");
273 return {};
274 }
275
276 // Cells
277 int const nCells = cells->GetNumberOfCells();
278 if(nCells > 0) {
279 if(!cells->IsStorage64Bit()) {
280 if(cells->CanConvertTo64BitStorage()) {
281 this->printWrn("Converting the cell array to 64-bit storage");
282 bool const success = cells->ConvertTo64BitStorage();
283 if(!success) {
284 this->printErr(
285 "Error converting the provided cell array to 64-bit storage");
286 return {};
287 }
288 } else {
289 this->printErr(
290 "Cannot convert the provided cell array to 64-bit storage");
291 return {};
292 }
293 }
294 auto connectivity = static_cast<vtkIdType *>(
295 ttkUtils::GetVoidPointer(cells->GetConnectivityArray()));
296 auto offsets = static_cast<vtkIdType *>(
297 ttkUtils::GetVoidPointer(cells->GetOffsetsArray()));
298
299 int status;
300 if(hasIndexArray) {
301 status
302 = triangulation->setStellarInputCells(nCells, connectivity, offsets);
303 } else {
304 status = triangulation->setInputCells(nCells, connectivity, offsets);
305 }
306
307 if(status != 0) {
308 this->printErr(
309 "Run the `vtkTetrahedralize` filter to resolve the issue.");
310 return {};
311 }
312 }
313
314 if(hasIndexArray) {
315 this->printMsg("Initializing Compact Triangulation", 1,
318 } else {
319 this->printMsg("Initializing Explicit Triangulation", 1,
322 }
323
324 return triangulation;
325}
326
328 ttkTriangulationFactory::CreateTriangulation(vtkDataSet *dataSet) {
329 switch(dataSet->GetDataObjectType()) {
330 case VTK_UNSTRUCTURED_GRID:
331 case VTK_POLY_DATA: {
332 return this->CreateExplicitTriangulation(
333 static_cast<vtkPointSet *>(dataSet));
334 }
335 case VTK_IMAGE_DATA: {
336 return this->CreateImplicitTriangulation((vtkImageData *)dataSet);
337 }
338 default: {
339 this->printErr("Unable to triangulate `"
340 + std::string(dataSet->GetClassName()) + "`");
341 }
342 }
343
344 return nullptr;
345}
346
348 int debugLevel, float cacheRatio, vtkDataSet *object) {
349 auto instance = &ttkTriangulationFactory::Instance;
350 instance->setDebugLevel(debugLevel);
351
352 auto key = ttkTriangulationFactory::GetKey(object);
353
354 ttk::Triangulation *triangulation{nullptr};
355 auto it = instance->registry.find(key);
356 if(it != instance->registry.end()) {
357 // object is the owner of the explicit or implicit triangulation
358 if(it->second.isValid(object)) {
359 instance->printMsg(
360 "Retrieving Existing Triangulation", ttk::debug::Priority::DETAIL);
361 triangulation = it->second.triangulation.get();
362 } else {
363 instance->printMsg(
364 "Existing Triangulation No Longer Valid", ttk::debug::Priority::DETAIL);
365 instance->registry.erase(key);
366 }
367 }
368
369 if(!triangulation && object->IsA("vtkImageData")) {
370 instance->FindImplicitTriangulation(
371 triangulation, static_cast<vtkImageData *>(object));
372 if(triangulation)
373 instance->printMsg("Retrieving Equivalent Implicit-Triangulation",
375 }
376
377 if(!triangulation) {
378 triangulation = instance->CreateTriangulation(object).release();
379 if(triangulation) {
380 instance->registry.emplace(std::piecewise_construct,
381 std::forward_as_tuple(key),
382 std::forward_as_tuple(object, triangulation));
383 }
384 }
385
386 instance->printMsg(
387 "# Registered Triangulations: " + std::to_string(instance->registry.size()),
389
390 if(triangulation) {
391 triangulation->setDebugLevel(debugLevel);
392 triangulation->setCacheSize(cacheRatio);
393 }
394
395 return triangulation;
396}
397
398int ttkTriangulationFactory::FindImplicitTriangulation(
399 ttk::Triangulation *&triangulation, vtkImageData *image) {
400
401 for(const auto &it : this->registry) {
402 if(it.second.owner->IsA("vtkImageData")) {
403 if(it.second.isValid(image)) {
404 triangulation = it.second.triangulation.get();
405 return 1;
406 }
407 }
408 }
409
410 return 0;
411}
412
414 switch(dataSet->GetDataObjectType()) {
415 case VTK_IMAGE_DATA: {
416 return (RegistryKey)dataSet;
417 }
418 default: {
419 auto cells = GetCells(dataSet);
420 if(cells)
421 return (RegistryKey)cells;
422 }
423 }
424
425 return 0;
426}
427
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
static RegistryKey GetKey(vtkDataSet *dataSet)
static ttk::Triangulation * GetTriangulation(int debugLevel, float cacheRatio, vtkDataSet *object)
static ttkTriangulationFactory Instance
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printMsg(const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
Definition Debug.h:118
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
double getElapsedTime()
Definition Timer.h:15
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
const char compactTriangulationIndex[]
Definition DataTypes.h:85
bool isValid(vtkDataSet *dataSet) const
RegistryValue(vtkDataSet *dataSet, ttk::Triangulation *triangulation_)
RegistryTriangulation triangulation
void Execute(vtkObject *, unsigned long ttkNotUsed(eventId), void *ttkNotUsed(callData)) override
static ttkOnDeleteCommand * New()
void Init(vtkDataSet *dataSet)
long long RegistryKey
std::unique_ptr< ttk::Triangulation > RegistryTriangulation
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)