TTK
Loading...
Searching...
No Matches
PersistenceDiagram.h
Go to the documentation of this file.
1
147
148#pragma once
149
150// base code includes
151#include <ApproximateTopology.h>
153#include <FTMTreePP.h>
156#include <ProgressiveTopology.h>
157#include <Triangulation.h>
158
159namespace ttk {
160
165 class PersistenceDiagram : virtual public Debug {
166
167 public:
168 enum class BACKEND {
169 FTM = 0,
174 };
175
177
178 inline void setBackend(const BACKEND be) {
179 this->BackEnd = be;
180 }
181
182 inline void setComputeMinSad(const bool data) {
183 this->dms_.setComputeMinSad(data);
184 }
185 inline void setComputeSadSad(const bool data) {
186 this->dms_.setComputeSadSad(data);
187 }
188 inline void setComputeSadMax(const bool data) {
189 this->dms_.setComputeSadMax(data);
190 }
191
196 template <typename scalarType, typename triangulationType>
197 void
198 augmentPersistenceDiagram(std::vector<PersistencePair> &persistencePairs,
199 const scalarType *const scalars,
200 const triangulationType *triangulation);
201
203 ftm::TreeType treeType,
204 const SimplexId vertexId) const;
205
206 void sortPersistenceDiagram(std::vector<PersistencePair> &diagram,
207 const SimplexId *const offsets) const;
208
209 template <typename scalarType>
211 ftm::FTMTreePP &tree,
212 const std::vector<
213 std::tuple<ttk::SimplexId, ttk::SimplexId, scalarType, bool>> &pairs,
214 std::vector<PersistencePair> &diagram) const;
215
224 template <typename scalarType, class triangulationType>
225 int execute(std::vector<PersistencePair> &CTDiagram,
226 const scalarType *inputScalars,
227 const size_t scalarsMTime,
228 const SimplexId *inputOffsets,
229 const triangulationType *triangulation,
230 const std::vector<bool> *updateMask = nullptr);
231
232 template <typename scalarType, class triangulationType>
233 int executeFTM(std::vector<PersistencePair> &CTDiagram,
234 const scalarType *inputScalars,
235 const SimplexId *inputOffsets,
236 const triangulationType *triangulation);
237
238 template <class triangulationType>
239 int executeProgressiveTopology(std::vector<PersistencePair> &CTDiagram,
240 const SimplexId *inputOffsets,
241 const triangulationType *triangulation);
242 template <typename scalarType, class triangulationType>
243 int executeApproximateTopology(std::vector<PersistencePair> &CTDiagram,
244 const scalarType *inputScalars,
245 const triangulationType *triangulation);
246
247 template <class triangulationType>
248 int executePersistentSimplex(std::vector<PersistencePair> &CTDiagram,
249 const SimplexId *inputOffsets,
250 const triangulationType *triangulation);
251
252 template <typename scalarType, class triangulationType>
253 int executeDiscreteMorseSandwich(std::vector<PersistencePair> &CTDiagram,
254 const scalarType *inputScalars,
255 const size_t scalarsMTime,
256 const SimplexId *inputOffsets,
257 const triangulationType *triangulation,
258 const std::vector<bool> *updateMask
259 = nullptr);
260
261 template <class triangulationType>
262 void checkProgressivityRequirement(const triangulationType *triangulation);
263
264 template <class triangulationType>
265 void checkManifold(const triangulationType *const triangulation);
266
267 inline void
290
291 inline void setOutputMonotonyOffsets(void *data) {
293 }
294 inline void setOutputOffsets(void *data) {
295 outputOffsets_ = data;
296 }
297 inline void setOutputScalars(void *data) {
298 outputScalars_ = data;
299 }
300 inline void setDeltaApproximate(double data) {
301 approxT_.setDelta(data);
302 }
303
304 protected:
305 bool IgnoreBoundary{false};
310
311 // int BackEnd{0};
313 // progressivity
316
319 bool IsResumable{false};
320 double TimeLimit{};
321
322 // Approximation
326 double Epsilon;
327 };
328} // namespace ttk
329
330template <typename scalarType>
332 ftm::FTMTreePP &tree,
333 const std::vector<
334 std::tuple<ttk::SimplexId, ttk::SimplexId, scalarType, bool>> &pairs,
335 std::vector<PersistencePair> &diagram) const {
336
337 const ttk::SimplexId numberOfPairs = pairs.size();
338 diagram.resize(numberOfPairs);
339 for(ttk::SimplexId i = 0; i < numberOfPairs; ++i) {
340 const ttk::SimplexId v0 = std::get<0>(pairs[i]);
341 const ttk::SimplexId v1 = std::get<1>(pairs[i]);
342 const bool type = std::get<3>(pairs[i]);
343
344 if(type == true) {
345 diagram[i] = PersistencePair{
347 v0, getNodeType(tree.getJoinTree(), ftm::TreeType::Join, v0), {}, {}},
349 v1, getNodeType(tree.getJoinTree(), ftm::TreeType::Join, v1), {}, {}},
350 0, true};
351 } else {
352 diagram[i] = PersistencePair{
354 v1,
356 {},
357 {}},
359 v0,
361 {},
362 {}},
363 2, true};
364 }
365 }
366
367 diagram.back().isFinite = false; // global extrema pair is infinite
368
369 return 0;
370}
371
372template <typename scalarType, typename triangulationType>
374 std::vector<PersistencePair> &persistencePairs,
375 const scalarType *const scalars,
376 const triangulationType *triangulation) {
377
378#ifdef TTK_ENABLE_OPENMP
379#pragma omp parallel for num_threads(threadNumber_)
380#endif // TTK_ENABLE_OPENMP
381 for(std::size_t i = 0; i < persistencePairs.size(); ++i) {
382 auto &pair{persistencePairs[i]};
383 triangulation->getVertexPoint(pair.birth.id, pair.birth.coords[0],
384 pair.birth.coords[1], pair.birth.coords[2]);
385 pair.birth.sfValue = scalars[pair.birth.id];
386 triangulation->getVertexPoint(pair.death.id, pair.death.coords[0],
387 pair.death.coords[1], pair.death.coords[2]);
388 pair.death.sfValue = scalars[pair.death.id];
389 }
390}
391
392template <typename scalarType, class triangulationType>
393int ttk::PersistenceDiagram::execute(std::vector<PersistencePair> &CTDiagram,
394 const scalarType *inputScalars,
395 const size_t scalarsMTime,
396 const SimplexId *inputOffsets,
397 const triangulationType *triangulation,
398 const std::vector<bool> *updateMask) {
399
401
402 checkProgressivityRequirement(triangulation);
403 checkManifold(triangulation);
404
405 Timer tm{};
406
407 switch(BackEnd) {
408 case BACKEND::PERSISTENT_SIMPLEX:
409 executePersistentSimplex(CTDiagram, inputOffsets, triangulation);
410 break;
411 case BACKEND::DISCRETE_MORSE_SANDWICH:
412 executeDiscreteMorseSandwich(CTDiagram, inputScalars, scalarsMTime,
413 inputOffsets, triangulation, updateMask);
414 break;
415 case BACKEND::PROGRESSIVE_TOPOLOGY:
416 executeProgressiveTopology(CTDiagram, inputOffsets, triangulation);
417 break;
418 case BACKEND::APPROXIMATE_TOPOLOGY:
419 executeApproximateTopology(CTDiagram, inputScalars, triangulation);
420 break;
421 case BACKEND::FTM:
422 executeFTM(CTDiagram, inputScalars, inputOffsets, triangulation);
423 break;
424 default:
425 printErr("No method was selected");
426 }
427
428 this->printMsg("Complete", 1.0, tm.getElapsedTime(), this->threadNumber_);
429
430 // augment persistence pairs with meta-data
431 augmentPersistenceDiagram(CTDiagram, inputScalars, triangulation);
432
433 // finally sort the diagram
434 sortPersistenceDiagram(CTDiagram, inputOffsets);
435
437
438 return 0;
439}
440
441template <class triangulationType>
443 std::vector<PersistencePair> &CTDiagram,
444 const SimplexId *inputOffsets,
445 const triangulationType *triangulation) {
446
447 Timer const tm{};
448 const auto dim = triangulation->getDimensionality();
449
450 std::vector<ttk::PersistentSimplexPairs::PersistencePair> pairs{};
451
452 psp_.setDebugLevel(this->debugLevel_);
453 psp_.setThreadNumber(this->threadNumber_);
454 psp_.computePersistencePairs(pairs, inputOffsets, *triangulation);
455 dms_.setInputOffsets(inputOffsets);
456
457 // convert PersistentSimplex pairs (with critical cells id) to PL
458 // pairs (with vertices id)
459
460#ifdef TTK_ENABLE_OPENMP
461#pragma omp parallel for num_threads(threadNumber_)
462#endif // TTK_ENABLE_OPENMP
463 for(size_t i = 0; i < pairs.size(); ++i) {
464 auto &pair{pairs[i]};
465 if(pair.type > 0) {
466 pair.birth = dms_.getCellGreaterVertex(
467 Cell{pair.type, pair.birth}, *triangulation);
468 }
469 if(pair.death != -1) {
470 pair.death = dms_.getCellGreaterVertex(
471 Cell{pair.type + 1, pair.death}, *triangulation);
472 }
473 }
474
475 CTDiagram.reserve(pairs.size() + 1);
476
477 // find the global maximum
478 const auto nVerts = triangulation->getNumberOfVertices();
479 const SimplexId globmax = std::distance(
480 inputOffsets, std::max_element(inputOffsets, inputOffsets + nVerts));
481
482 // convert pairs to the relevant format
483 for(const auto &p : pairs) {
484 const auto isFinite = (p.death >= 0);
485 const auto death = isFinite ? p.death : globmax;
486 if(p.type == 0) {
487 const auto deathType = (isFinite && dim > 1)
490 CTDiagram.emplace_back(PersistencePair{
492 CriticalVertex{death, deathType, {}, {}}, p.type, isFinite});
493 } else if(p.type == 1) {
494 const auto birthType
496 const auto deathType = (isFinite && dim == 3)
499 CTDiagram.emplace_back(PersistencePair{
500 CriticalVertex{p.birth, birthType, {}, {}},
501 CriticalVertex{death, deathType, {}, {}}, p.type, isFinite});
502 } else if(p.type == 2) {
503 CTDiagram.emplace_back(PersistencePair{
504 CriticalVertex{p.birth, CriticalType::Saddle2, {}, {}},
506 isFinite});
507 }
508 }
509
510 return 0;
511}
512
513template <typename scalarType, class triangulationType>
515 std::vector<PersistencePair> &CTDiagram,
516 const scalarType *inputScalars,
517 const size_t scalarsMTime,
518 const SimplexId *inputOffsets,
519 const triangulationType *triangulation,
520 const std::vector<bool> *updateMask) {
521
522 Timer const tm{};
523 const auto dim = triangulation->getDimensionality();
524
525 dms_.buildGradient(
526 inputScalars, scalarsMTime, inputOffsets, *triangulation, updateMask);
527 std::vector<DiscreteMorseSandwich::PersistencePair> dms_pairs{};
528 dms_.computePersistencePairs(
529 dms_pairs, inputOffsets, *triangulation, this->IgnoreBoundary);
530 CTDiagram.resize(dms_pairs.size());
531
532 // transform DiscreteMorseSandwich pairs (critical cells id) to PL
533 // pairs (vertices id)
534#ifdef TTK_ENABLE_OPENMP
535#pragma omp parallel for num_threads(threadNumber_)
536#endif // TTK_ENABLE_OPENMP
537 for(size_t i = 0; i < dms_pairs.size(); ++i) {
538 auto &pair{dms_pairs[i]};
539 if(pair.type > 0) {
540 pair.birth = dms_.getCellGreaterVertex(
541 Cell{pair.type, pair.birth}, *triangulation);
542 }
543 if(pair.death != -1) {
544 pair.death = dms_.getCellGreaterVertex(
545 Cell{pair.type + 1, pair.death}, *triangulation);
546 }
547 }
548
549 // find the global maximum
550 const auto nVerts = triangulation->getNumberOfVertices();
551 const SimplexId globmax = std::distance(
552 inputOffsets, std::max_element(inputOffsets, inputOffsets + nVerts));
553
554 // convert pairs to the relevant format
555#ifdef TTK_ENABLE_OPENMP
556#pragma omp parallel for num_threads(threadNumber_)
557#endif // TTK_ENABLE_OPENMP
558 for(size_t i = 0; i < dms_pairs.size(); ++i) {
559 const auto &p{dms_pairs[i]};
560 const auto isFinite = (p.death >= 0);
561 const auto death = isFinite ? p.death : globmax;
562
563 if(p.type == 0) {
564 const auto dtype = (isFinite && dim > 1) ? CriticalType::Saddle1
566 CTDiagram[i] = PersistencePair{
568 CriticalVertex{death, dtype, {}, {}}, p.type, isFinite};
569 } else if(p.type == 1) {
570 const auto btype
572 const auto dtype = (isFinite && dim == 3) ? CriticalType::Saddle2
574 CTDiagram[i] = PersistencePair{CriticalVertex{p.birth, btype, {}, {}},
575 CriticalVertex{death, dtype, {}, {}},
576 p.type, isFinite};
577 } else if(p.type == 2) {
578 const auto btype = (isFinite || dim == 3) ? CriticalType::Saddle2
580 CTDiagram[i] = PersistencePair{
581 CriticalVertex{p.birth, btype, {}, {}},
583 isFinite};
584 }
585 }
586
587 return 0;
588}
589
590template <typename scalarType, class triangulationType>
592 std::vector<PersistencePair> &CTDiagram,
593 const scalarType *inputScalars,
594 const triangulationType *triangulation) {
595
596 approxT_.setDebugLevel(debugLevel_);
597 approxT_.setThreadNumber(threadNumber_);
598 approxT_.setupTriangulation(const_cast<ttk::ImplicitTriangulation *>(
599 (const ImplicitTriangulation *)triangulation));
600 approxT_.setStartingResolutionLevel(StartingResolutionLevel);
601 approxT_.setStoppingResolutionLevel(StoppingResolutionLevel);
602 approxT_.setPreallocateMemory(true);
603 approxT_.setEpsilon(Epsilon);
604
605 std::vector<ApproximateTopology::PersistencePair> resultDiagram{};
606
607 approxT_.computeApproximatePD(
608 resultDiagram, inputScalars, (scalarType *)outputScalars_,
609 (SimplexId *)outputOffsets_, (int *)outputMonotonyOffsets_);
610
611 // create the final diagram
612 for(const auto &p : resultDiagram) {
613 if(p.pairType == 0) {
614 CTDiagram.emplace_back(PersistencePair{
616 CriticalVertex{p.death, CriticalType::Saddle1, {}, {}}, p.pairType,
617 true});
618 } else if(p.pairType == 2) {
619 CTDiagram.emplace_back(PersistencePair{
620 CriticalVertex{p.birth, CriticalType::Saddle2, {}, {}},
622 p.pairType, true});
623 } else if(p.pairType == -1) {
624 CTDiagram.emplace_back(PersistencePair{
627 p.pairType, false});
628 }
629 }
630
631 return 0;
632}
633
634template <class triangulationType>
636 std::vector<PersistencePair> &CTDiagram,
637 const SimplexId *inputOffsets,
638 const triangulationType *triangulation) {
639
640 progT_.setDebugLevel(debugLevel_);
641 progT_.setThreadNumber(threadNumber_);
642 progT_.setupTriangulation(const_cast<ttk::ImplicitTriangulation *>(
643 (const ImplicitTriangulation *)triangulation));
644 progT_.setStartingResolutionLevel(StartingResolutionLevel);
645 progT_.setStoppingResolutionLevel(StoppingResolutionLevel);
646 progT_.setTimeLimit(TimeLimit);
647 progT_.setIsResumable(IsResumable);
648 progT_.setPreallocateMemory(true);
649
650 std::vector<ProgressiveTopology::PersistencePair> resultDiagram{};
651
652 progT_.computeProgressivePD(resultDiagram, inputOffsets);
653
654 // create the final diagram
655 for(const auto &p : resultDiagram) {
656 if(p.pairType == 0) {
657 CTDiagram.emplace_back(PersistencePair{
659 CriticalVertex{p.death, CriticalType::Saddle1, {}, {}}, p.pairType,
660 true});
661 } else if(p.pairType == 2) {
662 CTDiagram.emplace_back(PersistencePair{
663 CriticalVertex{p.birth, CriticalType::Saddle2, {}, {}},
665 p.pairType, true});
666 } else if(p.pairType == -1) {
667 CTDiagram.emplace_back(PersistencePair{
669 CriticalVertex{p.death, CriticalType::Local_maximum, {}, {}}, 0,
670 false});
671 }
672 }
673
674 return 0;
675}
676
677template <typename scalarType, class triangulationType>
679 std::vector<PersistencePair> &CTDiagram,
680 const scalarType *inputScalars,
681 const SimplexId *inputOffsets,
682 const triangulationType *triangulation) {
683
684 contourTree_.setVertexScalars(inputScalars);
685 contourTree_.setTreeType(ftm::TreeType::Join_Split);
686 contourTree_.setVertexSoSoffsets(inputOffsets);
687 contourTree_.setSegmentation(false);
688 contourTree_.build<scalarType>(triangulation);
689
690 // get persistence pairs
691 std::vector<std::tuple<ttk::SimplexId, ttk::SimplexId, scalarType>> JTPairs;
692 std::vector<std::tuple<ttk::SimplexId, ttk::SimplexId, scalarType>> STPairs;
693 contourTree_.computePersistencePairs<scalarType>(JTPairs, true);
694 contourTree_.computePersistencePairs<scalarType>(STPairs, false);
695
696 // merge pairs
697 const auto JTSize = JTPairs.size();
698 const auto STSize = STPairs.size();
699 std::vector<std::tuple<ttk::SimplexId, ttk::SimplexId, scalarType, bool>>
700 CTPairs(JTSize + STSize);
701 for(size_t i = 0; i < JTSize; ++i) {
702 const auto &x = JTPairs[i];
703 CTPairs[i]
704 = std::make_tuple(std::get<0>(x), std::get<1>(x), std::get<2>(x), true);
705 }
706 for(size_t i = 0; i < STSize; ++i) {
707 const auto &x = STPairs[i];
708 CTPairs[JTSize + i]
709 = std::make_tuple(std::get<0>(x), std::get<1>(x), std::get<2>(x), false);
710 }
711
712 // remove the last pair which is present two times (global extrema pair)
713 if(!CTPairs.empty()) {
714 auto cmp =
715 [](
716 const std::tuple<ttk::SimplexId, ttk::SimplexId, scalarType, bool> &a,
717 const std::tuple<ttk::SimplexId, ttk::SimplexId, scalarType, bool> &b) {
718 return std::get<2>(a) < std::get<2>(b);
719 };
720
721 std::sort(CTPairs.begin(), CTPairs.end(), cmp);
722 CTPairs.erase(CTPairs.end() - 1);
723 }
724
725 // get persistence diagrams
726 computeCTPersistenceDiagram<scalarType>(contourTree_, CTPairs, CTDiagram);
727
728 return 0;
729}
730
731template <class triangulationType>
733 const triangulationType *ttkNotUsed(triangulation)) {
734
735 if((BackEnd == BACKEND::PROGRESSIVE_TOPOLOGY
736 || BackEnd == BACKEND::APPROXIMATE_TOPOLOGY)
737 && !std::is_same<ttk::ImplicitWithPreconditions, triangulationType>::value
738 && !std::is_same<ttk::ImplicitNoPreconditions, triangulationType>::value) {
739
740 printWrn("Explicit, Compact or Periodic triangulation detected.");
741 printWrn("Defaulting to the FTM backend.");
742
743 BackEnd = BACKEND::FTM;
744 }
745}
746
747template <class triangulationType>
749 const triangulationType *const triangulation) {
750
751 if(this->BackEnd != BACKEND::DISCRETE_MORSE_SANDWICH) {
752 return;
753 }
754
755 if(!triangulation->isManifold()) {
756 this->printWrn("Non-manifold data-set detected.");
757 this->printWrn("Defaulting to the Persistence Simplex backend.");
758
759 this->BackEnd = BACKEND::PERSISTENT_SIMPLEX;
760 }
761}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
AbstractTriangulation is an interface class that defines an interface for efficient traversal methods...
TTK processing package for progressive Topological Data Analysis.
virtual int setThreadNumber(const int threadNumber)
Definition BaseClass.h:80
Minimalist debugging class.
Definition Debug.h:88
int debugLevel_
Definition Debug.h:379
virtual int setDebugLevel(const int &debugLevel)
Definition Debug.cpp:147
TTK DiscreteMorseSandwich processing package.
void setComputeSadMax(const bool data)
void setComputeMinSad(const bool data)
void preconditionTriangulation(AbstractTriangulation *const data)
void setComputeSadSad(const bool data)
ImplicitTriangulation is a class that provides time and memory efficient traversal methods on triangu...
TTK processing package for the computation of persistence diagrams.
int executePersistentSimplex(std::vector< PersistencePair > &CTDiagram, const SimplexId *inputOffsets, const triangulationType *triangulation)
void setComputeSadSad(const bool data)
void setOutputMonotonyOffsets(void *data)
void sortPersistenceDiagram(std::vector< PersistencePair > &diagram, const SimplexId *const offsets) const
void preconditionTriangulation(AbstractTriangulation *triangulation)
ttk::ProgressiveTopology progT_
void setDeltaApproximate(double data)
void setBackend(const BACKEND be)
ttk::CriticalType getNodeType(ftm::FTMTree_MT *tree, ftm::TreeType treeType, const SimplexId vertexId) const
void checkProgressivityRequirement(const triangulationType *triangulation)
void setOutputScalars(void *data)
int executeFTM(std::vector< PersistencePair > &CTDiagram, const scalarType *inputScalars, const SimplexId *inputOffsets, const triangulationType *triangulation)
dcg::DiscreteGradient dcg_
int executeProgressiveTopology(std::vector< PersistencePair > &CTDiagram, const SimplexId *inputOffsets, const triangulationType *triangulation)
int execute(std::vector< PersistencePair > &CTDiagram, const scalarType *inputScalars, const size_t scalarsMTime, const SimplexId *inputOffsets, const triangulationType *triangulation, const std::vector< bool > *updateMask=nullptr)
ttk::ApproximateTopology approxT_
DiscreteMorseSandwich dms_
int executeDiscreteMorseSandwich(std::vector< PersistencePair > &CTDiagram, const scalarType *inputScalars, const size_t scalarsMTime, const SimplexId *inputOffsets, const triangulationType *triangulation, const std::vector< bool > *updateMask=nullptr)
void setComputeMinSad(const bool data)
void setOutputOffsets(void *data)
void augmentPersistenceDiagram(std::vector< PersistencePair > &persistencePairs, const scalarType *const scalars, const triangulationType *triangulation)
Complete a ttk::DiagramType instance with scalar field values (useful for persistence) and 3D coordin...
void setComputeSadMax(const bool data)
int executeApproximateTopology(std::vector< PersistencePair > &CTDiagram, const scalarType *inputScalars, const triangulationType *triangulation)
PersistentSimplexPairs psp_
void checkManifold(const triangulationType *const triangulation)
int computeCTPersistenceDiagram(ftm::FTMTreePP &tree, const std::vector< std::tuple< ttk::SimplexId, ttk::SimplexId, scalarType, bool > > &pairs, std::vector< PersistencePair > &diagram) const
Textbook algorithm to find persistence pairs.
void preconditionTriangulation(AbstractTriangulation *const data)
Preprocess all the required connectivity requests on the triangulation.
TTK processing package for progressive Topological Data Analysis.
TTK discreteGradient processing package.
int setThreadNumber(const int n) override
Definition FTMTree_CT.h:86
FTMTree_MT * getSplitTree()
Definition FTMTree_CT.h:51
int setDebugLevel(const int &d) override
Definition FTMTree_CT.h:79
FTMTree_MT * getJoinTree()
Definition FTMTree_CT.h:47
void preconditionTriangulation(AbstractTriangulation *tri, const bool preproc=true)
Definition FTMTree_CT.h:72
The Topology ToolKit.
CriticalType
default value for critical index
Definition DataTypes.h:80
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)