TTK
Loading...
Searching...
No Matches
ttkPersistenceDiagramClustering.cpp
Go to the documentation of this file.
1#include <ttkMacros.h>
4#include <ttkUtils.h>
5
6#include <vtkCellData.h>
7#include <vtkDataArray.h>
8#include <vtkDoubleArray.h>
9#include <vtkFieldData.h>
10#include <vtkFloatArray.h>
11#include <vtkInformation.h>
12#include <vtkIntArray.h>
13#include <vtkMultiBlockDataSet.h>
14#include <vtkNew.h>
15#include <vtkPointData.h>
16#include <vtkUnsignedCharArray.h>
17#include <vtkUnstructuredGrid.h>
18
19using namespace ttk;
20
22
24 SetNumberOfInputPorts(1);
25 SetNumberOfOutputPorts(3);
26}
27
29 int port, vtkInformation *info) {
30 if(port == 0) {
31 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
32 } else {
33 return 0;
34 }
35 return 1;
36}
37
39 int port, vtkInformation *info) {
40 if(port == 0 || port == 1 || port == 2) {
41 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
42 } else {
43 return 0;
44 }
45 return 1;
46}
47
49 needUpdate_ = true;
50 ttkAlgorithm::Modified();
51}
52
54 vtkInformation *ttkNotUsed(request),
55 vtkInformationVector **inputVector,
56 vtkInformationVector *outputVector) {
57
58 auto blocks = vtkMultiBlockDataSet::GetData(inputVector[0], 0);
59
60 // Flat storage for diagrams extracted from blocks
61 std::vector<vtkUnstructuredGrid *> input;
62
63 // Number of input diagrams
64 int numInputs = 0;
65
66 if(blocks != nullptr) {
67 numInputs = blocks->GetNumberOfBlocks();
68 input.resize(numInputs);
69 for(int i = 0; i < numInputs; ++i) {
70 input[i] = vtkUnstructuredGrid::SafeDownCast(blocks->GetBlock(i));
71 if(this->GetMTime() < input[i]->GetMTime()) {
72 needUpdate_ = true;
73 }
74 }
75 }
76
77 if(numInputs == 0) {
78 this->printErr("No input detected");
79 return 0;
80 }
81 if(numInputs != 2 and NonMatchingWeight != 1.0) {
82 // We need to modify the update of the barycenter by taking into account the
83 // non-matching weight
84 printWrn("Custom weight only supported for pairwise distance.");
85 printWrn("Non-matching weight is set to 1.");
87 }
88
89 // Get output pointers
90 auto output_clusters = vtkMultiBlockDataSet::GetData(outputVector, 0);
91 auto output_centroids = vtkMultiBlockDataSet::GetData(outputVector, 1);
92 auto output_matchings = vtkMultiBlockDataSet::GetData(outputVector, 2);
93
94 if(needUpdate_) {
95 // clear data before computation
96 intermediateDiagrams_ = {};
97 all_matchings_ = {};
98 final_centroids_ = {};
99
100 intermediateDiagrams_.resize(numInputs);
101 all_matchings_.resize(3);
102
103 // store the persistence of every min-max global pair
104 std::vector<double> max_persistences(numInputs);
105
106 for(int i = 0; i < numInputs; i++) {
107 auto &diag{this->intermediateDiagrams_[i]};
108 const auto ret = VTUToDiagram(diag, input[i], *this);
109 if(ret < 0) {
110 this->printErr("Could not read Persistence Diagram");
111 return 0;
112 }
113 if(this->NumberOfClusters > 1) {
114 // duplicate the global min-max pair in 2: one min-saddle pair and
115 // one saddle-max pair
116 diag[0].death.type = ttk::CriticalType::Saddle1;
117 // store the saddle max pair at the vector end
118 diag.emplace_back(diag[0]);
119 diag.back().birth.type = ttk::CriticalType::Saddle1;
120 diag.back().death.type = ttk::CriticalType::Local_maximum;
121 }
122 max_persistences[i] = diag[0].persistence();
123 }
124
125 this->max_dimension_total_
126 = *std::max_element(max_persistences.begin(), max_persistences.end());
127
128 if(this->Method == METHOD::PROGRESSIVE) {
129
130 if(!UseInterruptible) {
131 TimeLimit = 999999999;
132 }
133
134 inv_clustering_ = this->execute(
135 intermediateDiagrams_, final_centroids_, all_matchings_);
136 needUpdate_ = false;
137
138 } else if(this->Method == METHOD::AUCTION) {
139
140 final_centroids_.resize(1);
141 inv_clustering_.resize(numInputs);
142 for(int i_input = 0; i_input < numInputs; i_input++) {
143 inv_clustering_[i_input] = 0;
144 }
145 PersistenceDiagramBarycenter pdBarycenter{};
146
147 const auto wassersteinMetric = std::to_string(WassersteinMetric);
148 pdBarycenter.setWasserstein(wassersteinMetric);
149 pdBarycenter.setMethod(2);
150 pdBarycenter.setNumberOfInputs(numInputs);
151 pdBarycenter.setDeterministic(Deterministic);
152 pdBarycenter.setUseProgressive(UseProgressive);
153 pdBarycenter.setDebugLevel(debugLevel_);
154 pdBarycenter.setThreadNumber(threadNumber_);
155 pdBarycenter.setAlpha(Alpha);
156 pdBarycenter.setLambda(Lambda);
157 pdBarycenter.setNonMatchingWeight(NonMatchingWeight);
158 pdBarycenter.setDeltaLim(DeltaLim);
159 pdBarycenter.execute(
160 intermediateDiagrams_, final_centroids_[0], all_matchings_);
161
162 needUpdate_ = false;
163 }
164 }
165
166 outputClusteredDiagrams(output_clusters, input, this->intermediateDiagrams_,
167 this->all_matchings_, this->inv_clustering_,
168 this->DisplayMethod, this->Spacing,
169 this->max_dimension_total_);
170 outputCentroids(output_centroids, this->final_centroids_,
171 this->all_matchings_, input[0], this->DisplayMethod,
172 this->Spacing, this->max_dimension_total_);
174 output_matchings, this->NumberOfClusters, this->intermediateDiagrams_,
175 this->all_matchings_, this->final_centroids_, this->inv_clustering_,
176 this->DisplayMethod, this->Spacing, this->max_dimension_total_);
177
178 return 1;
179}
180
181static void addCostsAsFieldData(vtkUnstructuredGrid *vtu,
182 const double minSadCost,
183 const double sadSadCost,
184 const double sadMaxCost) {
185
186 // add global matchings cost as FieldData (only 1 tuple per array)
187 vtkNew<vtkDoubleArray> minSad{};
188 minSad->SetName("MinSaddleCost");
189 minSad->SetNumberOfTuples(1);
190 minSad->SetTuple1(0, minSadCost);
191 vtu->GetFieldData()->AddArray(minSad);
192
193 vtkNew<vtkDoubleArray> sadSad{};
194 sadSad->SetName("SaddleSaddleCost");
195 sadSad->SetNumberOfTuples(1);
196 sadSad->SetTuple1(0, sadSadCost);
197 vtu->GetFieldData()->AddArray(sadSad);
198
199 vtkNew<vtkDoubleArray> sadMax{};
200 sadMax->SetName("SaddleMaxCost");
201 sadMax->SetNumberOfTuples(1);
202 sadMax->SetTuple1(0, sadMaxCost);
203 vtu->GetFieldData()->AddArray(sadMax);
204
205 vtkNew<vtkDoubleArray> wass{};
206 wass->SetName("WassersteinDistance");
207 wass->SetNumberOfTuples(1);
208 wass->SetTuple1(0, minSadCost + sadSadCost + sadMaxCost);
209 vtu->GetFieldData()->AddArray(wass);
210}
211
213 vtkMultiBlockDataSet *output,
214 const std::vector<vtkUnstructuredGrid *> &diagsVTU,
215 const std::vector<ttk::DiagramType> &diags,
216 const std::vector<std::vector<std::vector<ttk::MatchingType>>>
217 &matchingsPerCluster,
218 const std::vector<int> &inv_clustering,
219 const DISPLAY dm,
220 const double spacing,
221 const double max_persistence) const {
222
223 // index of diagram in its cluster
224 std::vector<int> diagIdInClust{};
225 // number of diagrams per cluster
226 std::vector<int> clustSize{};
227
228 // prep work for displaying diagrams as cluster stars
229 if(dm == DISPLAY::STARS) {
230 // total number of clusters
231 const auto nClusters
232 = 1 + *std::max_element(inv_clustering.begin(), inv_clustering.end());
233 clustSize.resize(nClusters, 0);
234 diagIdInClust.resize(diagsVTU.size());
235 for(size_t i = 0; i < inv_clustering.size(); ++i) {
236 auto &diagsInClust = clustSize[inv_clustering[i]];
237 diagIdInClust[i] = diagsInClust;
238 diagsInClust++;
239 }
240 }
241
242 output->SetNumberOfBlocks(diagsVTU.size());
243
244 for(size_t i = 0; i < diagsVTU.size(); ++i) {
245 vtkNew<vtkUnstructuredGrid> vtu{};
246 vtu->ShallowCopy(diagsVTU[i]);
247
248 // put clusterId on FieldData (only 1 tuple)
249 vtkNew<vtkIntArray> clusterId{};
250 clusterId->SetName("ClusterId");
251 clusterId->SetNumberOfComponents(1);
252 clusterId->SetNumberOfTuples(1);
253 clusterId->Fill(inv_clustering[i]);
254 vtu->GetFieldData()->AddArray(clusterId);
255
256 // add ClusterID also on PointData?
257 vtkNew<vtkIntArray> cidPD{};
258 cidPD->SetName("ClusterID");
259 cidPD->SetNumberOfComponents(1);
260 cidPD->SetNumberOfTuples(vtu->GetNumberOfPoints());
261 cidPD->Fill(inv_clustering[i]);
262 vtu->GetPointData()->AddArray(cidPD);
263
264 // add Persistence also on PointData?
265 vtkNew<vtkDoubleArray> pointPers{};
266 pointPers->SetName("Persistence");
267 pointPers->SetNumberOfTuples(vtu->GetNumberOfPoints());
268 vtu->GetPointData()->AddArray(pointPers);
269
270 // diagonal uses two existing points
271 const auto persArray = vtu->GetCellData()->GetArray("Persistence");
272 for(int j = 0; j < vtu->GetNumberOfCells() - 1; ++j) {
273 const auto pers = persArray->GetTuple1(j);
274 pointPers->SetTuple1(2 * j + 0, pers);
275 pointPers->SetTuple1(2 * j + 1, pers);
276 }
277 pointPers->SetTuple1(vtu->GetNumberOfPoints() - 1,
278 persArray->GetTuple1(vtu->GetNumberOfCells() - 1));
279
280 const auto cid = inv_clustering[i];
281 const auto &matchings{matchingsPerCluster[cid][i]};
282 double minSadCost{}, sadSadCost{}, sadMaxCost{};
283 const auto &diag{diags[i]};
284
285 for(size_t j = 0; j < matchings.size(); ++j) {
286 const auto &m{matchings[j]};
287 const auto bidderId{std::get<0>(m)};
288
289 // avoid out-of-bound accesses
290 if(bidderId >= static_cast<ttk::SimplexId>(diag.size())) {
291 this->printWrn("Out-of-bounds access averted");
292 continue;
293 }
294
295 if(bidderId > 0) {
296 const auto &p1{diag[bidderId]};
297 if(p1.birth.type == ttk::CriticalType::Local_minimum) {
298 minSadCost += std::get<2>(m);
299 } else if(p1.birth.type == ttk::CriticalType::Saddle1
300 && p1.death.type == ttk::CriticalType::Saddle2) {
301 sadSadCost += std::get<2>(m);
302 } else if(p1.death.type == ttk::CriticalType::Local_maximum) {
303 sadMaxCost += std::get<2>(m);
304 }
305 }
306 }
307
308 addCostsAsFieldData(vtu, minSadCost, sadSadCost, sadMaxCost);
309
310 // translate diagram back to canonical representation
311 ResetDiagramPosition(vtu, *this);
312
313 if(dm == DISPLAY::MATCHINGS && spacing > 0) {
314 // translate diagrams along the Z axis
315 const std::array<double, 3> trans{0, 0, i == 0 ? -spacing : spacing};
316 TranslateDiagram(vtu, trans);
317 output->SetBlock(i, vtu);
318 } else if(dm == DISPLAY::STARS && spacing > 0) {
319 const auto c = inv_clustering[i];
320 const auto angle = 2.0 * M_PI * static_cast<double>(diagIdInClust[i])
321 / static_cast<double>(clustSize[c]);
322 // translate diagrams in the XY plane
323 const std::array<double, 3> trans{
324 3.0 * (spacing + 0.2) * max_persistence * c
325 + spacing * max_persistence * std::cos(angle) + 0.2,
326 spacing * max_persistence * std::sin(angle), 0};
327 TranslateDiagram(vtu, trans);
328 output->SetBlock(i, vtu);
329
330 } else {
331 // add diagram to output multi-block dataset
332 output->SetBlock(i, vtu);
333 }
334 }
335}
336
338 vtkMultiBlockDataSet *output,
339 const std::vector<DiagramType> &final_centroids,
340 const std::vector<std::vector<std::vector<ttk::MatchingType>>>
341 &matchingsPerCluster,
342 vtkUnstructuredGrid *const someInputDiag,
343 const DISPLAY dm,
344 const double spacing,
345 const double max_persistence) const {
346
347 if(final_centroids.size() != matchingsPerCluster.size()) {
348 this->printWrn("Inconsistent matchings vector size");
349 }
350
351 const auto da{
352 someInputDiag->GetCellData()->GetArray(ttk::PersistenceBirthName)};
353 const auto dim{static_cast<int>(someInputDiag->GetCellData()
355 ->GetRange()[1])
356 + 1};
357
358 for(size_t i = 0; i < final_centroids.size(); ++i) {
359 vtkNew<vtkUnstructuredGrid> vtu{};
360 DiagramToVTU(vtu, final_centroids[i], da, *this, dim, false);
361
362 // put clusterId on FieldData (only 1 tuple)
363 vtkNew<vtkIntArray> cid{};
364 cid->SetName("ClusterId");
365 cid->SetNumberOfTuples(1);
366 cid->SetTuple1(0, i);
367 vtu->GetFieldData()->AddArray(cid);
368
369 // add ClusterID also on PointData?
370 vtkNew<vtkIntArray> cidPD{};
371 cidPD->SetName("ClusterID");
372 cidPD->SetNumberOfComponents(1);
373 cidPD->SetNumberOfTuples(vtu->GetNumberOfPoints());
374 cidPD->Fill(i);
375 vtu->GetPointData()->AddArray(cidPD);
376
377 // add Persistence also on PointData?
378 vtkNew<vtkDoubleArray> pointPers{};
379 pointPers->SetName("Persistence");
380 pointPers->SetNumberOfTuples(vtu->GetNumberOfPoints());
381 vtu->GetPointData()->AddArray(pointPers);
382
383 // diagonal uses two existing points
384 for(int j = 0; j < vtu->GetNumberOfCells() - 1; ++j) {
385 const auto persArray = vtu->GetCellData()->GetArray("Persistence");
386 const auto pers = persArray->GetTuple1(j);
387 pointPers->SetTuple1(2 * j + 0, pers);
388 pointPers->SetTuple1(2 * j + 1, pers);
389 }
390
391 double minSadCost{}, sadSadCost{}, sadMaxCost{};
392
393 for(const auto &matchingsPerDiag : matchingsPerCluster[i]) {
394 for(const auto &m : matchingsPerDiag) {
395 const auto goodId{std::get<1>(m)};
396 const auto &p0{final_centroids[i][goodId]};
397 if(p0.birth.type == ttk::CriticalType::Local_minimum) {
398 minSadCost += std::get<2>(m);
399 } else if(p0.birth.type == ttk::CriticalType::Saddle1
400 && p0.death.type == ttk::CriticalType::Saddle2) {
401 sadSadCost += std::get<2>(m);
402 } else if(p0.death.type == ttk::CriticalType::Local_maximum) {
403 sadMaxCost += std::get<2>(m);
404 }
405 }
406 }
407
408 addCostsAsFieldData(vtu, minSadCost, sadSadCost, sadMaxCost);
409
410 if(dm == DISPLAY::STARS && spacing > 0) {
411 // shift centroid along the X axis
412 const std::array<double, 3> trans{
413 3.0 * (spacing + 0.2) * max_persistence * i, 0, 0};
414 TranslateDiagram(vtu, trans);
415 output->SetBlock(i, vtu);
416
417 } else {
418 // add centroid to output multi-block dataset
419 output->SetBlock(i, vtu);
420 }
421 }
422}
423
425 vtkMultiBlockDataSet *output,
426 const size_t nClusters,
427 const std::vector<DiagramType> &diags,
428 const std::vector<std::vector<std::vector<ttk::MatchingType>>>
429 &matchingsPerCluster,
430 const std::vector<DiagramType> &centroids,
431 const std::vector<int> &inv_clustering,
432 const DISPLAY dm,
433 const double spacing,
434 const double max_persistence) const {
435
436 // index of diagram in its cluster
437 std::vector<int> diagIdInClust{};
438 // number of diagrams per cluster
439 std::vector<int> clustSize{};
440
441 // prep work for displaying diagrams as cluster stars
442 if(dm == DISPLAY::STARS) {
443 clustSize.resize(nClusters, 0);
444 diagIdInClust.resize(diags.size());
445 for(size_t i = 0; i < inv_clustering.size(); ++i) {
446 auto &diagsInClust = clustSize[inv_clustering[i]];
447 diagIdInClust[i] = diagsInClust;
448 diagsInClust++;
449 }
450 }
451
452 // count the number of bidders per centroid pair
453 // (when with only 1 cluster and 2 diagrams)
454 std::vector<int> matchings_count(centroids[0].size());
455 std::vector<int> count_to_good{};
456
457 for(size_t i = 0; i < diags.size(); ++i) {
458 const auto cid = inv_clustering[i];
459 const auto &diag{diags[i]};
460 const auto &matchings{matchingsPerCluster[cid][i]};
461
462 vtkNew<vtkUnstructuredGrid> matchingsGrid{};
463
464 const auto nCells{matchings.size()};
465 const auto nPoints{2 * matchings.size()};
466
467 vtkNew<vtkPoints> points{};
468 points->SetNumberOfPoints(nPoints);
469 matchingsGrid->SetPoints(points);
470
471 // point data
472 vtkNew<vtkIntArray> diagIdVerts{};
473 diagIdVerts->SetName("DiagramID");
474 diagIdVerts->SetNumberOfTuples(nPoints);
475 matchingsGrid->GetPointData()->AddArray(diagIdVerts);
476
477 vtkNew<vtkIntArray> pointId{};
478 pointId->SetName("PointID");
479 pointId->SetNumberOfTuples(nPoints);
480 matchingsGrid->GetPointData()->AddArray(pointId);
481
482 // cell data
483 vtkNew<vtkIntArray> diagIdCells{};
484 diagIdCells->SetName("DiagramID");
485 diagIdCells->SetNumberOfTuples(nCells);
486 matchingsGrid->GetCellData()->AddArray(diagIdCells);
487
488 vtkNew<vtkIntArray> clusterId{};
489 clusterId->SetName("ClusterID");
490 clusterId->SetNumberOfTuples(nCells);
491 clusterId->Fill(cid);
492 matchingsGrid->GetCellData()->AddArray(clusterId);
493
494 vtkNew<vtkDoubleArray> matchCost{};
495 matchCost->SetName("Cost");
496 matchCost->SetNumberOfTuples(nCells);
497 matchingsGrid->GetCellData()->AddArray(matchCost);
498
499 vtkNew<vtkIntArray> pairType{};
500 pairType->SetName("PairType");
501 pairType->SetNumberOfTuples(nCells);
502 matchingsGrid->GetCellData()->AddArray(pairType);
503
504 vtkNew<vtkIntArray> isDiagonal{};
505 isDiagonal->SetName("IsDiagonal");
506 isDiagonal->SetNumberOfTuples(nCells);
507 matchingsGrid->GetCellData()->AddArray(isDiagonal);
508
509 for(size_t j = 0; j < matchings.size(); ++j) {
510 const auto &m{matchings[j]};
511 const auto bidderId{std::get<0>(m)};
512 const auto goodId{std::get<1>(m)};
513
514 // avoid out-of-bound accesses
515 if(goodId >= static_cast<ttk::SimplexId>(centroids[cid].size())
516 || bidderId >= static_cast<ttk::SimplexId>(diag.size())) {
517 this->printWrn("Out-of-bounds access averted");
518 continue;
519 }
520
521 if(nClusters == 1) {
522 matchings_count[goodId] += 1;
523 count_to_good.push_back(goodId);
524 }
525
526 const auto &p0{centroids[cid][goodId]};
527 std::array<double, 3> coords1{};
528 std::array<double, 3> coords0{p0.birth.sfValue, p0.death.sfValue, 0};
529
530 if(bidderId >= 0) {
531 const auto &p{diag[bidderId]};
532 coords1[0] = p.birth.sfValue;
533 coords1[1] = p.death.sfValue;
534 coords1[2] = 0;
535 isDiagonal->SetTuple1(j, 0);
536 } else {
537 const double diagonal_projection
538 = (p0.birth.sfValue + p0.death.sfValue) / 2.0;
539 coords1[0] = diagonal_projection;
540 coords1[1] = diagonal_projection;
541 coords1[2] = 0;
542 isDiagonal->SetTuple1(j, 1);
543 }
544
545 if(dm == DISPLAY::STARS && spacing > 0) {
546 const auto angle = 2.0 * M_PI * static_cast<double>(diagIdInClust[i])
547 / static_cast<double>(clustSize[cid]);
548 const auto shift
549 = 3.0 * (std::abs(spacing) + 0.2) * max_persistence * cid;
550 coords0[0] += shift;
551 coords1[0] += shift + spacing * max_persistence * std::cos(angle);
552 coords1[1] += spacing * max_persistence * std::sin(angle);
553
554 } else if(dm == DISPLAY::MATCHINGS) {
555 coords1[2] = (diags.size() == 2 && i == 0) ? -spacing : spacing;
556 }
557
558 points->SetPoint(2 * j + 0, coords0.data());
559 points->SetPoint(2 * j + 1, coords1.data());
560 std::array<vtkIdType, 2> ids{
561 2 * static_cast<vtkIdType>(j) + 0,
562 2 * static_cast<vtkIdType>(j) + 1,
563 };
564 matchingsGrid->InsertNextCell(VTK_LINE, 2, ids.data());
565
566 diagIdCells->SetTuple1(j, i);
567 matchCost->SetTuple1(j, std::get<2>(m));
568 diagIdVerts->SetTuple1(2 * j + 0, i);
569 diagIdVerts->SetTuple1(2 * j + 1, i);
570 pointId->SetTuple1(2 * j + 0, goodId);
571 pointId->SetTuple1(2 * j + 1, bidderId);
572
573 if(bidderId >= 0) {
574 const auto &p1{diag[bidderId]};
575 pairType->SetTuple1(j, p1.isFinite ? p1.dim : -1);
576 } else {
577 pairType->SetTuple1(j, p0.dim);
578 }
579 }
580
581 addCostsAsFieldData(matchingsGrid, this->distances[0], this->distances[1],
582 this->distances[2]);
583
584 // add diagram matchings to multi-block
585 output->SetBlock(i, matchingsGrid);
586 }
587
588 // add matchings number
589 if(nClusters == 1 && diags.size() == 2) {
590 size_t nPrevCells{};
591 for(size_t i = 0; i < diags.size(); ++i) {
592 vtkNew<vtkIntArray> matchNumber{};
593 matchNumber->SetName("MatchNumber");
594 const auto matchings
595 = vtkUnstructuredGrid::SafeDownCast(output->GetBlock(i));
596 const auto nCells = matchings->GetNumberOfCells();
597 matchNumber->SetNumberOfTuples(nCells);
598 for(int j = 0; j < nCells; j++) {
599 const auto goodId = count_to_good[j + (i == 0 ? 0 : nPrevCells)];
600 matchNumber->SetTuple1(j, matchings_count[goodId]);
601 }
602 matchings->GetCellData()->AddArray(matchNumber);
603 nPrevCells = nCells;
604 }
605 }
606}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
#define M_PI
Definition Os.h:50
TTK processing package for the computation of Wasserstein barycenters and K-Means clusterings of a se...
int FillInputPortInformation(int port, vtkInformation *info) override
void outputCentroids(vtkMultiBlockDataSet *output, const std::vector< ttk::DiagramType > &final_centroids, const std::vector< std::vector< std::vector< ttk::MatchingType > > > &matchingsPerCluster, vtkUnstructuredGrid *const someInputDiag, const DISPLAY dm, const double spacing, const double max_persistence) const
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillOutputPortInformation(int port, vtkInformation *info) override
void outputClusteredDiagrams(vtkMultiBlockDataSet *output, const std::vector< vtkUnstructuredGrid * > &diagsVTU, const std::vector< ttk::DiagramType > &diags, const std::vector< std::vector< std::vector< ttk::MatchingType > > > &matchingsPerCluster, const std::vector< int > &inv_clustering, const DISPLAY dm, const double spacing, const double max_persistence) const
void outputMatchings(vtkMultiBlockDataSet *output, const size_t nClusters, const std::vector< ttk::DiagramType > &diags, const std::vector< std::vector< std::vector< ttk::MatchingType > > > &matchingsPerCluster, const std::vector< ttk::DiagramType > &centroids, const std::vector< int > &inv_clustering, const ttkPersistenceDiagramClustering::DISPLAY dm, const double spacing, const double max_persistence) const
int debugLevel_
Definition Debug.h:379
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
std::vector< int > execute(std::vector< DiagramType > &intermediateDiagrams, std::vector< DiagramType > &centroids, std::vector< std::vector< std::vector< MatchingType > > > &all_matchings)
The Topology ToolKit.
const char PersistencePairTypeName[]
Definition DataTypes.h:73
const char PersistenceBirthName[]
Definition DataTypes.h:68
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkPersistenceDiagramClustering)
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 ResetDiagramPosition(vtkUnstructuredGrid *const diagram, const ttk::Debug &dbg)
Translate back a canonical diagram into its original position.