TTK
Loading...
Searching...
No Matches
ttkPersistenceDiagramUtils.cpp
Go to the documentation of this file.
1#include <ttkMacros.h>
3#include <ttkUtils.h>
4
5#include <vtkCellData.h>
6#include <vtkDoubleArray.h>
7#include <vtkFloatArray.h>
8#include <vtkPointData.h>
9#include <vtkThreshold.h>
10#include <vtkTransform.h>
11#include <vtkTransformFilter.h>
12#include <vtkUnsignedCharArray.h>
13#include <vtkUnstructuredGrid.h>
14#include <vtkVersionMacros.h> // for VTK_VERSION_CHECK
15
17 vtkUnstructuredGrid *vtu,
18 const ttk::Debug &dbg) {
19
20 const auto pd = vtu->GetPointData();
21 const auto cd = vtu->GetCellData();
22 const auto points = vtu->GetPoints();
23
24 if(pd == nullptr) {
25 dbg.printWrn("VTU diagram with NULL Point Data");
26 return -1;
27 }
28 if(cd == nullptr) {
29 dbg.printWrn("VTU diagram with NULL Cell Data");
30 return -2;
31 }
32 if(points == nullptr) {
33 dbg.printWrn("VTU with no points");
34 return -3;
35 }
36
37 // cell data
38 const auto pairId = ttkSimplexIdTypeArray::SafeDownCast(
40 const auto pairType
41 = vtkIntArray::SafeDownCast(cd->GetArray(ttk::PersistencePairTypeName));
42 const auto pairPers = cd->GetArray(ttk::PersistenceName);
43 const auto birthScalars = cd->GetArray(ttk::PersistenceBirthName);
44 const auto isFinite = cd->GetArray(ttk::PersistenceIsFinite);
45
46 // point data
47 const auto vertexId = ttkSimplexIdTypeArray::SafeDownCast(
48 pd->GetArray(ttk::VertexScalarFieldName));
49 const auto critType
50 = vtkIntArray::SafeDownCast(pd->GetArray(ttk::PersistenceCriticalTypeName));
51 const auto coords = vtkFloatArray::SafeDownCast(
53
54 const bool embed = coords == nullptr;
55
56 if(pairId == nullptr) {
57 dbg.printWrn("Missing PairIdentifier cell data array");
58 return -5;
59 }
60 if(pairType == nullptr) {
61 dbg.printWrn("Missing PairType cell data array");
62 return -6;
63 }
64 if(pairPers == nullptr) {
65 dbg.printWrn("Missing Persistence cell data array");
66 return -7;
67 }
68 if(vertexId == nullptr) {
69 dbg.printWrn("Missing ttkVertexScalarField point data array");
70 return -8;
71 }
72 if(critType == nullptr) {
73 dbg.printWrn("Missing CriticalType point data array");
74 return -9;
75 }
76 if(birthScalars == nullptr) {
77 dbg.printWrn("Missing Birth cell data array");
78 return -10;
79 }
80 if(isFinite == nullptr) {
81 dbg.printWrn("Missing IsFinite cell data array");
82 return -12;
83 }
84
85 int nPairs = pairId->GetNumberOfTuples();
86
87 // compact pairIds in [0, nPairs - 1] (diagonal excepted)
88 for(int i = 0; i < nPairs; i++) {
89 if(pairId->GetTuple1(i) != -1) {
90 pairId->SetTuple1(i, i);
91 } else {
92 // detect diagram diagonal
93 nPairs -= 1;
94 }
95 }
96
97 if(nPairs < 1) {
98 dbg.printWrn("Diagram has no pairs");
99 return -4;
100 }
101
102 diagram.resize(nPairs);
103
104 // skip diagonal cell (corresponding points already dealt with)
105 for(int i = 0; i < nPairs; ++i) {
106
107 const auto pId = pairId->GetValue(i);
108 // skip diagram diagonal if present
109 if(pId == -1) {
110 continue;
111 }
112
113 const auto i0 = 2 * i + 0;
114 const auto i1 = 2 * i + 1;
115
116 const auto v0 = vertexId->GetValue(i0);
117 const auto v1 = vertexId->GetValue(i1);
118 const auto ct0 = static_cast<ttk::CriticalType>(critType->GetValue(i0));
119 const auto ct1 = static_cast<ttk::CriticalType>(critType->GetValue(i1));
120
121 const auto pType = pairType->GetValue(i);
122 const auto pers = pairPers->GetTuple1(i);
123 const auto birth = birthScalars->GetTuple1(i);
124 const auto isFin = static_cast<bool>(isFinite->GetTuple1(i));
125
126 std::array<float, 3> coordsBirth{}, coordsDeath{};
127 // no vtkPoints::GetPoint() taking a float array, have to do the
128 // conversion by hand...
129 std::array<double, 3> tmp{};
130
131 if(embed) {
132 points->GetPoint(i0, tmp.data());
133 coordsBirth = {static_cast<float>(tmp[0]), static_cast<float>(tmp[1]),
134 static_cast<float>(tmp[2])};
135 points->GetPoint(i1, tmp.data());
136 coordsDeath = {static_cast<float>(tmp[0]), static_cast<float>(tmp[1]),
137 static_cast<float>(tmp[2])};
138 } else {
139 coords->GetTuple(i0, tmp.data());
140 coordsBirth = {static_cast<float>(tmp[0]), static_cast<float>(tmp[1]),
141 static_cast<float>(tmp[2])};
142 coords->GetTuple(i1, tmp.data());
143 coordsDeath = {static_cast<float>(tmp[0]), static_cast<float>(tmp[1]),
144 static_cast<float>(tmp[2])};
145 }
146
147 // put pairs in diagram
148 diagram[i] = ttk::PersistencePair{
149 ttk::CriticalVertex{v0, ct0, birth, coordsBirth},
150 ttk::CriticalVertex{v1, ct1, birth + pers, coordsDeath}, pType, isFin};
151 }
152
153 return 0;
154}
155
156int DiagramToVTU(vtkUnstructuredGrid *vtu,
157 const ttk::DiagramType &diagram,
158 vtkDataArray *const inputScalars,
159 const ttk::Debug &dbg,
160 const int dim,
161 const bool embedInDomain) {
162
163 if(diagram.empty()) {
164 dbg.printWrn("Empty diagram");
165 return -1;
166 }
167
168 const auto pd = vtu->GetPointData();
169 const auto cd = vtu->GetCellData();
170
171 if(pd == nullptr || cd == nullptr) {
172 dbg.printWrn("Grid has no point data or no cell data");
173 return -2;
174 }
175
176 // point data arrays
177
178 vtkNew<ttkSimplexIdTypeArray> vertsId{};
179 vertsId->SetName(ttk::VertexScalarFieldName);
180 vertsId->SetNumberOfTuples(2 * diagram.size());
181 pd->AddArray(vertsId);
182
183 vtkNew<vtkIntArray> critType{};
184 critType->SetName(ttk::PersistenceCriticalTypeName);
185 critType->SetNumberOfTuples(2 * diagram.size());
186 pd->AddArray(critType);
187
188 vtkNew<vtkFloatArray> coordsScalars{};
189
190 if(!embedInDomain) {
191 coordsScalars->SetNumberOfComponents(3);
192 coordsScalars->SetName(ttk::PersistenceCoordinatesName);
193 coordsScalars->SetNumberOfTuples(2 * diagram.size());
194 pd->AddArray(coordsScalars);
195 }
196
197 // cell data arrays
198
199 vtkNew<ttkSimplexIdTypeArray> pairsId{};
200 pairsId->SetName(ttk::PersistencePairIdentifierName);
201 pairsId->SetNumberOfTuples(diagram.size());
202 cd->AddArray(pairsId);
203
204 vtkNew<vtkIntArray> pairsDim{};
205 pairsDim->SetName(ttk::PersistencePairTypeName);
206 pairsDim->SetNumberOfTuples(diagram.size());
207 cd->AddArray(pairsDim);
208
209 vtkSmartPointer<vtkDataArray> const persistence{inputScalars->NewInstance()};
210 persistence->SetName(ttk::PersistenceName);
211 persistence->SetNumberOfTuples(diagram.size());
212 cd->AddArray(persistence);
213
214 vtkSmartPointer<vtkDataArray> const birthScalars{inputScalars->NewInstance()};
215 birthScalars->SetName(ttk::PersistenceBirthName);
216 birthScalars->SetNumberOfTuples(diagram.size());
217 cd->AddArray(birthScalars);
218
219 vtkNew<vtkUnsignedCharArray> isFinite{};
220 isFinite->SetName(ttk::PersistenceIsFinite);
221 isFinite->SetNumberOfTuples(diagram.size());
222 cd->AddArray(isFinite);
223
224 // grid
225
226 vtkNew<vtkPoints> points{};
227 points->SetNumberOfPoints(2 * diagram.size());
228 vtkNew<vtkIdTypeArray> offsets{}, connectivity{};
229 offsets->SetNumberOfComponents(1);
230 offsets->SetNumberOfTuples(diagram.size() + 1);
231 connectivity->SetNumberOfComponents(1);
232 connectivity->SetNumberOfTuples(2 * diagram.size());
233
234#ifdef TTK_ENABLE_OPENMP
235#pragma omp parallel for num_threads(dbg.getThreadNumber())
236#endif // TTK_ENABLE_OPENMP
237 for(size_t i = 0; i < diagram.size(); ++i) {
238 const auto &pair{diagram[i]};
239 const auto i0{2 * i + 0}, i1{2 * i + 1};
240 if(embedInDomain) {
241 points->SetPoint(
242 i0, pair.birth.coords[0], pair.birth.coords[1], pair.birth.coords[2]);
243 points->SetPoint(
244 i1, pair.death.coords[0], pair.death.coords[1], pair.death.coords[2]);
245 } else {
246 points->SetPoint(i0, pair.birth.sfValue, pair.birth.sfValue, 0);
247 points->SetPoint(i1, pair.birth.sfValue, pair.death.sfValue, 0);
248 }
249
250 connectivity->SetTuple1(i0, i0);
251 connectivity->SetTuple1(i1, i1);
252 offsets->SetTuple1(i, 2 * i);
253
254 // point data
255 vertsId->SetTuple1(i0, pair.birth.id);
256 vertsId->SetTuple1(i1, pair.death.id);
257 critType->SetTuple1(i0, static_cast<ttk::SimplexId>(pair.birth.type));
258 critType->SetTuple1(i1, static_cast<ttk::SimplexId>(pair.death.type));
259
260 if(!embedInDomain) {
261 coordsScalars->SetTuple3(
262 i0, pair.birth.coords[0], pair.birth.coords[1], pair.birth.coords[2]);
263 coordsScalars->SetTuple3(
264 i1, pair.death.coords[0], pair.death.coords[1], pair.death.coords[2]);
265 }
266
267 // cell data
268 pairsId->SetTuple1(i, i);
269 persistence->SetTuple1(i, pair.persistence());
270 birthScalars->SetTuple1(i, pair.birth.sfValue);
271 isFinite->SetTuple1(i, pair.isFinite);
272 pairsDim->SetTuple1(
273 i, (pair.dim == 2 && pair.isFinite) ? dim - 1 : pair.dim);
274 }
275 offsets->SetTuple1(diagram.size(), connectivity->GetNumberOfTuples());
276
277 vtkNew<vtkCellArray> cells{};
278 cells->SetData(offsets, connectivity);
279 vtu->SetPoints(points);
280 vtu->SetCells(VTK_LINE, cells);
281
282 if(!embedInDomain) {
283 // highest birth (last pair)
284 const auto lastPair = std::max_element(diagram.begin(), diagram.end());
285 // add diagonal (first point -> last birth/penultimate point)
286 std::array<vtkIdType, 2> diag{
287 0, 2 * std::distance(diagram.begin(), lastPair)};
288 vtu->InsertNextCell(VTK_LINE, 2, diag.data());
289 pairsId->InsertTuple1(diagram.size(), -1);
290 pairsDim->InsertTuple1(diagram.size(), -1);
291 isFinite->InsertTuple1(diagram.size(), false);
292 // persistence of global min-max pair
293 const auto maxPersistence = diagram[0].persistence();
294 persistence->InsertTuple1(diagram.size(), 2 * maxPersistence);
295 // birth == death == 0
296 birthScalars->InsertTuple1(diagram.size(), 0);
297 }
298
299 return 0;
300}
301
302int ProjectDiagramInsideDomain(vtkUnstructuredGrid *const inputDiagram,
303 vtkUnstructuredGrid *const outputDiagram,
304 const ttk::Debug &dbg) {
305
306 ttk::Timer tm{};
307
308 // use vtkThreshold to remove diagonal (PairIdentifier == -1)
309 vtkNew<vtkThreshold> threshold{};
310 threshold->SetInputDataObject(0, inputDiagram);
311 threshold->SetInputArrayToProcess(0, 0, 0,
312 vtkDataObject::FIELD_ASSOCIATION_CELLS,
314#if VTK_VERSION_NUMBER < VTK_VERSION_CHECK(9, 2, 0)
315 threshold->ThresholdByUpper(0);
316#else
317 threshold->SetThresholdFunction(vtkThreshold::THRESHOLD_UPPER);
318 threshold->SetUpperThreshold(0);
319#endif
320
321 threshold->Update();
322
323 auto diagonalLess = threshold->GetOutput();
324 auto diagonalLessData = diagonalLess->GetPointData();
325
326 const auto critCoordinates = vtkFloatArray::SafeDownCast(
327 diagonalLessData->GetAbstractArray(ttk::PersistenceCoordinatesName));
328
329 // set new points from Coordinates array
330 vtkNew<vtkFloatArray> coords{};
331 coords->DeepCopy(critCoordinates);
332 coords->SetName("Points");
333 diagonalLess->GetPoints()->SetData(coords);
334 diagonalLessData->RemoveArray(ttk::PersistenceCoordinatesName);
335
336 outputDiagram->ShallowCopy(diagonalLess);
337
338 // don't forget to forward the Field Data
339 outputDiagram->GetFieldData()->ShallowCopy(inputDiagram->GetFieldData());
340
341 dbg.printMsg("Projected Persistence Diagram inside domain", 1.0,
342 tm.getElapsedTime(), dbg.getThreadNumber());
343
344 return 0;
345}
346
347template <typename dataType>
348void getCoords(vtkPoints *points,
349 const dataType *const births,
350 const dataType *const perss,
351 vtkIdType nPoints,
352 const int nThreads) {
353
354#ifdef TTK_ENABLE_OPENMP
355#pragma omp parallel for num_threads(nThreads)
356#endif // TTK_ENABLE_OPENMP
357 for(int i = 0; i < nPoints / 2; ++i) {
358 std::array<float, 3> pt0{
359 static_cast<float>(births[i]),
360 static_cast<float>(births[i]),
361 0,
362 };
363 std::array<float, 3> pt1{
364 static_cast<float>(births[i]),
365 static_cast<float>(births[i] + perss[i]),
366 0,
367 };
368 points->SetPoint(2 * i + 0, pt0.data());
369 points->SetPoint(2 * i + 1, pt1.data());
370 }
371
372 TTK_FORCE_USE(nThreads);
373}
374
375int ProjectDiagramIn2D(vtkUnstructuredGrid *const inputDiagram,
376 vtkUnstructuredGrid *const outputDiagram,
377 const ttk::Debug &dbg) {
378
379 ttk::Timer tm{};
380
381 outputDiagram->ShallowCopy(inputDiagram);
382
383 auto pointData = outputDiagram->GetPointData();
384
385 auto birth = inputDiagram->GetCellData()->GetArray(ttk::PersistenceBirthName);
386 auto pers = inputDiagram->GetCellData()->GetArray(ttk::PersistenceName);
387
388 if(birth == nullptr || pers == nullptr) {
389 dbg.printWrn("Missing Birth or Persistence arrays");
390 return 1;
391 }
392
393 // generate a new `Coordinates` pointData array
394 vtkNew<vtkFloatArray> coords{};
395 coords->DeepCopy(inputDiagram->GetPoints()->GetData());
396 coords->SetName(ttk::PersistenceCoordinatesName);
397 pointData->AddArray(coords);
398
399 vtkNew<vtkPoints> points{};
400 const auto nPoints = inputDiagram->GetNumberOfPoints();
401 points->SetNumberOfPoints(nPoints);
402
403 if(birth->GetNumberOfTuples() != nPoints / 2
404 || pers->GetNumberOfTuples() != nPoints / 2) {
405 dbg.printWrn("Wrong number of tuples for Birth or Persistence arrays");
406 return 2;
407 }
408
409 switch(birth->GetDataType()) {
410 vtkTemplateMacro(getCoords(points, ttkUtils::GetPointer<VTK_TT>(birth),
411 ttkUtils::GetPointer<VTK_TT>(pers), nPoints,
412 dbg.getThreadNumber()));
413 }
414
415 outputDiagram->SetPoints(points);
416
417 // add diagonal(first point -> last birth/penultimate point)
418 std::array<vtkIdType, 2> diag{0, 2 * (outputDiagram->GetNumberOfCells() - 1)};
419 outputDiagram->InsertNextCell(VTK_LINE, 2, diag.data());
420
421 // add diagonal data
422 auto cellData = outputDiagram->GetCellData();
423 auto pairIdentifierScalars = vtkIntArray::SafeDownCast(
424 cellData->GetArray(ttk::PersistencePairIdentifierName));
425 auto extremumIndexScalars = vtkIntArray::SafeDownCast(
426 cellData->GetArray(ttk::PersistencePairTypeName));
427 auto persistenceScalars = cellData->GetArray(ttk::PersistenceName);
428 auto birthScalars = cellData->GetArray(ttk::PersistenceBirthName);
429 auto isFinite = cellData->GetArray(ttk::PersistenceIsFinite);
430
431 pairIdentifierScalars->InsertNextTuple1(-1);
432 extremumIndexScalars->InsertNextTuple1(-1);
433 isFinite->InsertNextTuple1(0);
434 // 2 * persistence of min-max pair
435 persistenceScalars->InsertNextTuple1(2 * persistenceScalars->GetTuple1(0));
436 // birth == death == 0
437 birthScalars->InsertNextTuple1(0);
438
439 // don't forget to forward the Field Data
440 outputDiagram->GetFieldData()->ShallowCopy(inputDiagram->GetFieldData());
441
442 dbg.printMsg("Projected Persistence Diagram back to 2D", 1.0,
443 tm.getElapsedTime(), dbg.getThreadNumber());
444
445 return 0;
446}
447
448int TranslateDiagram(vtkUnstructuredGrid *const diagram,
449 const std::array<double, 3> &trans) {
450
451 vtkNew<vtkUnstructuredGrid> tmp{};
452 tmp->ShallowCopy(diagram);
453
454 vtkNew<vtkTransform> tr{};
455 tr->Translate(trans.data());
456
457 vtkNew<vtkTransformFilter> trf{};
458 trf->SetTransform(tr);
459 trf->SetInputData(tmp);
460 trf->Update();
461
462 diagram->ShallowCopy(trf->GetOutputDataObject(0));
463
464 return 0;
465}
466
467int ResetDiagramPosition(vtkUnstructuredGrid *const diagram,
468 const ttk::Debug &dbg) {
469
470 const bool embedded
471 = diagram->GetPointData()->GetArray(ttk::PersistenceCoordinatesName)
472 == nullptr;
473
474 if(embedded) {
475 dbg.printWrn("Cannot reset embedded diagram position");
476 return 1;
477 }
478
479 // position of first point in diagram
480 std::array<double, 3> pos{};
481 diagram->GetPoint(diagram->GetCell(0)->GetPointId(0), pos.data());
482
483 // birth value of the first cell in the diagram
484 const auto firstBirth{
485 diagram->GetCellData()->GetArray(ttk::PersistenceBirthName)->GetTuple1(0)};
486
487 if((pos[0] != pos[1] && pos[0] != firstBirth) || pos[2] != 0) {
488 vtkNew<vtkUnstructuredGrid> tmp{};
489 tmp->ShallowCopy(diagram);
490
491 vtkNew<vtkTransform> tr{};
492 tr->Translate(firstBirth - pos[0], firstBirth - pos[1], -pos[2]);
493
494 vtkNew<vtkTransformFilter> trf{};
495 trf->SetTransform(tr);
496 trf->SetInputData(tmp);
497 trf->Update();
498
499 diagram->ShallowCopy(trf->GetOutputDataObject(0));
500
501 dbg.printMsg("Diagram reset to initial position");
502 }
503
504 return 0;
505}
#define TTK_FORCE_USE(x)
Force the compiler to use the function/method parameter.
Definition BaseClass.h:57
int getThreadNumber() const
Definition BaseClass.h:76
Minimalist debugging class.
Definition Debug.h:88
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
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
const char PersistenceName[]
Definition DataTypes.h:72
CriticalType
default value for critical index
Definition DataTypes.h:80
const char PersistenceCoordinatesName[]
Definition DataTypes.h:70
const char PersistencePairTypeName[]
Definition DataTypes.h:73
std::vector< PersistencePair > DiagramType
Persistence Diagram type as a vector of Persistence pairs.
const char PersistenceCriticalTypeName[]
Definition DataTypes.h:67
const char PersistenceIsFinite[]
Definition DataTypes.h:74
const char VertexScalarFieldName[]
default name for vertex scalar field
Definition DataTypes.h:35
const char PersistencePairIdentifierName[]
Definition DataTypes.h:71
const char PersistenceBirthName[]
Definition DataTypes.h:68
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
int ProjectDiagramInsideDomain(vtkUnstructuredGrid *const inputDiagram, vtkUnstructuredGrid *const outputDiagram, const ttk::Debug &dbg)
Generate the spatial embedding of a given Persistence Diagram.
int DiagramToVTU(vtkUnstructuredGrid *vtu, const ttk::DiagramType &diagram, vtkDataArray *const inputScalars, const ttk::Debug &dbg, const int dim, const bool embedInDomain)
Converts a Persistence Diagram in the ttk::DiagramType format to the VTK Unstructured Grid format (as...
int TranslateDiagram(vtkUnstructuredGrid *const diagram, const std::array< double, 3 > &trans)
Translate a diagram to a new position.
int VTUToDiagram(ttk::DiagramType &diagram, vtkUnstructuredGrid *vtu, const ttk::Debug &dbg)
Converts a Persistence Diagram in the VTK Unstructured Grid format (as generated by the ttkPersistenc...
int ProjectDiagramIn2D(vtkUnstructuredGrid *const inputDiagram, vtkUnstructuredGrid *const outputDiagram, const ttk::Debug &dbg)
Generate the 2D embedding of a given Persistence Diagram.
void getCoords(vtkPoints *points, const dataType *const births, const dataType *const perss, vtkIdType nPoints, const int nThreads)
int ResetDiagramPosition(vtkUnstructuredGrid *const diagram, const ttk::Debug &dbg)
Translate back a canonical diagram into its original position.