TTK
Loading...
Searching...
No Matches
ScalarFieldCriticalPoints.h
Go to the documentation of this file.
1
43
44#pragma once
45
46#include <map>
47
48// base code includes
49#include <ProgressiveTopology.h>
50#include <Triangulation.h>
51#include <UnionFind.h>
52
53namespace ttk {
54
55 class ScalarFieldCriticalPoints : public virtual Debug {
56
57 public:
59
60 enum class BACKEND { GENERIC = 0, PROGRESSIVE_TOPOLOGY = 1 };
61
75 template <class triangulationType = AbstractTriangulation>
76 int execute(const SimplexId *const offsets,
77 const triangulationType *triangulation);
78
79 template <class triangulationType = AbstractTriangulation>
80 int executeLegacy(const SimplexId *const offsets,
81 const triangulationType *triangulation);
82
83 template <class triangulationType = AbstractTriangulation>
84 int executeProgressive(const SimplexId *const offsets,
85 const triangulationType *triangulation);
86
87 template <class triangulationType>
88 void checkProgressivityRequirement(const triangulationType *triangulation);
89
90 template <class triangulationType = AbstractTriangulation>
92 const SimplexId vertexId,
93 const SimplexId *const offsets,
94 const triangulationType *triangulation,
95 bool &isLowerOnBoundary,
96 bool &isUpperOnBoundary,
97 std::vector<std::vector<ttk::SimplexId>> *upperComponents,
98 std::vector<std::vector<ttk::SimplexId>> *lowerComponents) const;
99
100 template <class triangulationType = AbstractTriangulation>
101 char
102 getCriticalType(const SimplexId &vertexId,
103 const SimplexId *const offsets,
104 const triangulationType *triangulation,
105 std::vector<std::vector<ttk::SimplexId>> *upperComponents
106 = nullptr,
107 std::vector<std::vector<ttk::SimplexId>> *lowerComponents
108 = nullptr) const;
109
110 char getCriticalType(const SimplexId &vertexId,
111 const SimplexId *const offsets,
112 const std::vector<std::pair<SimplexId, SimplexId>>
113 &vertexLinkEdgeList) const;
114
115 inline void setDomainDimension(const int &dimension) {
116 dimension_ = dimension;
117 }
118
119 inline void
120 setOutput(std::vector<std::pair<SimplexId, char>> *criticalPoints) {
121 criticalPoints_ = criticalPoints;
122 }
123
124 inline void
126 // pre-condition functions
127 if(triangulation) {
128 triangulation->preconditionVertexNeighbors();
129 triangulation->preconditionVertexStars();
130 triangulation->preconditionBoundaryVertices();
131 }
132
133 setDomainDimension(triangulation->getDimensionality());
134 setVertexNumber(triangulation->getNumberOfVertices());
135 }
136
138 const std::vector<std::vector<std::pair<SimplexId, SimplexId>>>
139 *edgeList) {
140 vertexLinkEdgeLists_ = edgeList;
141 }
142
146 int setVertexNumber(const SimplexId &vertexNumber) {
147 vertexNumber_ = vertexNumber;
148 return 0;
149 }
150
151 void setNonManifold(const bool b) {
153 }
154
155 void displayStats();
156
157 protected:
160 const std::vector<std::vector<std::pair<SimplexId, SimplexId>>>
162 std::vector<std::pair<SimplexId, char>> *criticalPoints_{};
163
165
166 // progressive
171 bool IsResumable{false};
172 double TimeLimit{};
173 };
174} // namespace ttk
175
176// template functions
177template <class triangulationType>
179 const SimplexId *const offsets, const triangulationType *triangulation) {
180
181 checkProgressivityRequirement(triangulation);
182
183 switch(BackEnd) {
184
186 this->executeProgressive(offsets, triangulation);
187 break;
188
189 case BACKEND::GENERIC:
190 this->executeLegacy(offsets, triangulation);
191 break;
192
193 default:
194 printErr("No method was selected");
195 }
196
198 return 0;
199}
200
201template <class triangulationType>
203 const SimplexId *const offsets, const triangulationType *triangulation) {
204
205 progT_.setDebugLevel(debugLevel_);
206 progT_.setThreadNumber(threadNumber_);
207 progT_.setupTriangulation(const_cast<ttk::ImplicitTriangulation *>(
208 (const ImplicitTriangulation *)triangulation));
209 progT_.setStartingResolutionLevel(StartingResolutionLevel);
210 progT_.setStoppingResolutionLevel(StoppingResolutionLevel);
211 progT_.setTimeLimit(TimeLimit);
212 progT_.setIsResumable(IsResumable);
213 progT_.setPreallocateMemory(true);
214
215 progT_.computeProgressiveCP(criticalPoints_, offsets);
216
217 displayStats();
218 return 0;
219}
220
221template <class triangulationType>
223 const SimplexId *const offsets, const triangulationType *triangulation) {
224
225 // check the consistency of the variables -- to adapt
226#ifndef TTK_ENABLE_KAMIKAZE
227 if((!dimension_) && ((!triangulation) || (triangulation->isEmpty())))
228 return -1;
229 if((!vertexNumber_) && ((!triangulation) || (triangulation->isEmpty())))
230 return -2;
231 if((!vertexLinkEdgeLists_) && (!triangulation))
232 return -4;
233 if(!criticalPoints_)
234 return -5;
235#endif
236
237 if(triangulation) {
238 vertexNumber_ = triangulation->getNumberOfVertices();
239 dimension_ = triangulation->getCellVertexNumber(0) - 1;
240 }
241
242 printMsg("Extracting critical points...");
243
244 Timer t;
245
246 std::vector<char> vertexTypes(vertexNumber_, (char)(CriticalType::Regular));
247
248#ifdef TTK_ENABLE_OPENMP
249 int const chunkSize
250 = std::max(1000, (int)vertexNumber_ / (threadNumber_ * 100));
251#endif
252
253 if(triangulation) {
254
255#ifdef TTK_ENABLE_OPENMP
256#pragma omp parallel for schedule(dynamic, chunkSize) num_threads(threadNumber_)
257#endif
258 for(SimplexId i = 0; i < (SimplexId)vertexNumber_; i++) {
259#ifdef TTK_ENABLE_MPI
260 if(!isRunningWithMPI()
261 || (isRunningWithMPI()
262 && (triangulation->getVertexRank(i) == ttk::MPIrank_)))
263#endif // TTK_ENABLE_MPI
264 {
265 vertexTypes[i] = getCriticalType(i, offsets, triangulation);
266 }
267 }
268 } else if(vertexLinkEdgeLists_) {
269 // legacy implementation
270#ifdef TTK_ENABLE_OPENMP
271#pragma omp parallel for schedule(dynamic, chunkSize) num_threads(threadNumber_)
272#endif
273 for(SimplexId i = 0; i < (SimplexId)vertexNumber_; i++) {
274 // can't use MPI there since triangulation (~> rankArray) is nullptr
275 vertexTypes[i] = getCriticalType(i, offsets, (*vertexLinkEdgeLists_)[i]);
276 }
277 }
278
279 SimplexId minimumNumber = 0, maximumNumber = 0, saddleNumber = 0,
280 oneSaddleNumber = 0, twoSaddleNumber = 0, monkeySaddleNumber = 0;
281
282 // debug msg
283 if(debugLevel_ >= (int)debug::Priority::INFO) {
284 if(dimension_ == 3) {
285 for(SimplexId i = 0; i < vertexNumber_; i++) {
286 switch(vertexTypes[i]) {
287
288 case(char)(CriticalType::Local_minimum):
289 minimumNumber++;
290 break;
291
292 case(char)(CriticalType::Saddle1):
293 oneSaddleNumber++;
294 break;
295
296 case(char)(CriticalType::Saddle2):
297 twoSaddleNumber++;
298 break;
299
300 case(char)(CriticalType::Local_maximum):
301 maximumNumber++;
302 break;
303
304 case(char)(CriticalType::Degenerate):
305 monkeySaddleNumber++;
306 break;
307 }
308 }
309 } else if(dimension_ == 2) {
310 for(SimplexId i = 0; i < vertexNumber_; i++) {
311 switch(vertexTypes[i]) {
312
313 case(char)(CriticalType::Local_minimum):
314 minimumNumber++;
315 break;
316
317 case(char)(CriticalType::Saddle1):
318 saddleNumber++;
319 break;
320
321 case(char)(CriticalType::Local_maximum):
322 maximumNumber++;
323 break;
324
325 case(char)(CriticalType::Degenerate):
326 monkeySaddleNumber++;
327 break;
328 }
329 }
330 }
331
332 {
333 std::vector<std::vector<std::string>> stats;
334 stats.push_back({" #Minima", std::to_string(minimumNumber)});
335 if(dimension_ == 3) {
336 stats.push_back({" #1-saddles", std::to_string(oneSaddleNumber)});
337 stats.push_back({" #2-saddles", std::to_string(twoSaddleNumber)});
338 }
339 if(dimension_ == 2) {
340 stats.push_back({" #Saddles", std::to_string(saddleNumber)});
341 }
342 stats.push_back({" #Multi-saddles", std::to_string(monkeySaddleNumber)});
343 stats.push_back({" #Maxima", std::to_string(maximumNumber)});
344
345 printMsg(stats);
346 }
347 }
348
349 // prepare the output
350 criticalPoints_->clear();
351 criticalPoints_->reserve(vertexNumber_);
352 for(SimplexId i = 0; i < vertexNumber_; i++) {
353 if(vertexTypes[i] != (char)(CriticalType::Regular)) {
354 criticalPoints_->emplace_back(i, vertexTypes[i]);
355 }
356 }
357
358 printMsg("Processed " + std::to_string(vertexNumber_) + " vertices", 1,
359 t.getElapsedTime(), threadNumber_);
360
361 return 0;
362}
363
364template <class triangulationType>
366 const SimplexId vertexId,
367 const SimplexId *const offsets,
368 const triangulationType *triangulation,
369 bool &isLowerOnBoundary,
370 bool &isUpperOnBoundary,
371 std::vector<std::vector<ttk::SimplexId>> *upperComponents,
372 std::vector<std::vector<ttk::SimplexId>> *lowerComponents) const {
373
374 SimplexId const neighborNumber
375 = triangulation->getVertexNeighborNumber(vertexId);
376 std::vector<SimplexId> lowerNeighbors, upperNeighbors;
377
378 for(SimplexId i = 0; i < neighborNumber; i++) {
379 SimplexId neighborId = 0;
380 triangulation->getVertexNeighbor(vertexId, i, neighborId);
381
382 if(offsets[neighborId] < offsets[vertexId]) {
383 lowerNeighbors.push_back(neighborId);
384 if(dimension_ == 3) {
385 if(triangulation->isVertexOnBoundary(neighborId))
386 isLowerOnBoundary = true;
387 }
388 }
389
390 // upper link
391 if(offsets[neighborId] > offsets[vertexId]) {
392 upperNeighbors.push_back(neighborId);
393 if(dimension_ == 3) {
394 if(triangulation->isVertexOnBoundary(neighborId))
395 isUpperOnBoundary = true;
396 }
397 }
398 }
399
400 // now do the actual work
401 std::vector<UnionFind> lowerSeeds(lowerNeighbors.size());
402 std::vector<UnionFind *> lowerList(lowerNeighbors.size());
403 std::vector<UnionFind> upperSeeds(upperNeighbors.size());
404 std::vector<UnionFind *> upperList(upperNeighbors.size());
405
406 for(SimplexId i = 0; i < (SimplexId)lowerSeeds.size(); i++) {
407 lowerList[i] = &(lowerSeeds[i]);
408 }
409 for(SimplexId i = 0; i < (SimplexId)upperSeeds.size(); i++) {
410 upperList[i] = &(upperSeeds[i]);
411 }
412
413 SimplexId const vertexStarSize = triangulation->getVertexStarNumber(vertexId);
414
415 for(SimplexId i = 0; i < vertexStarSize; i++) {
416 SimplexId cellId = 0;
417 triangulation->getVertexStar(vertexId, i, cellId);
418
419 SimplexId const cellSize = triangulation->getCellVertexNumber(cellId);
420 for(SimplexId j = 0; j < cellSize; j++) {
421 SimplexId neighborId0 = -1;
422 triangulation->getCellVertex(cellId, j, neighborId0);
423
424 if(neighborId0 != vertexId) {
425 // we are on the link
426
427 bool const lower0 = offsets[neighborId0] < offsets[vertexId];
428
429 // connect it to everybody except himself and vertexId
430 for(SimplexId k = j + 1; k < cellSize; k++) {
431
432 SimplexId neighborId1 = -1;
433 triangulation->getCellVertex(cellId, k, neighborId1);
434
435 if((neighborId1 != neighborId0) && (neighborId1 != vertexId)) {
436
437 bool const lower1 = offsets[neighborId1] < offsets[vertexId];
438
439 std::vector<SimplexId> *neighbors = &lowerNeighbors;
440 std::vector<UnionFind *> *seeds = &lowerList;
441
442 if(!lower0) {
443 neighbors = &upperNeighbors;
444 seeds = &upperList;
445 }
446
447 if(lower0 == lower1) {
448 // connect their union-find sets!
449 SimplexId lowerId0 = -1, lowerId1 = -1;
450 for(SimplexId l = 0; l < (SimplexId)neighbors->size(); l++) {
451 if((*neighbors)[l] == neighborId0) {
452 lowerId0 = l;
453 }
454 if((*neighbors)[l] == neighborId1) {
455 lowerId1 = l;
456 }
457 }
458 if((lowerId0 != -1) && (lowerId1 != -1)) {
459 (*seeds)[lowerId0] = UnionFind::makeUnion(
460 (*seeds)[lowerId0], (*seeds)[lowerId1]);
461 (*seeds)[lowerId1] = (*seeds)[lowerId0];
462 }
463 }
464 }
465 }
466 }
467 }
468 }
469
470 // let's remove duplicates now
471
472 // update the UF if necessary
473 for(SimplexId i = 0; i < (SimplexId)lowerList.size(); i++)
474 lowerList[i] = lowerList[i]->find();
475 for(SimplexId i = 0; i < (SimplexId)upperList.size(); i++)
476 upperList[i] = upperList[i]->find();
477
478 std::unordered_map<UnionFind *, std::vector<ttk::SimplexId>>::iterator it;
479 std::unordered_map<UnionFind *, std::vector<ttk::SimplexId>>
480 upperComponentId{};
481 std::unordered_map<UnionFind *, std::vector<ttk::SimplexId>>
482 lowerComponentId{};
483
484 // We retrieve the lower and upper components if we want them
485 for(ttk::SimplexId i = 0; i < (SimplexId)upperNeighbors.size(); i++) {
486 it = upperComponentId.find(upperList[i]);
487 if(it != upperComponentId.end()) {
488 upperComponentId[upperList[i]].push_back(upperNeighbors[i]);
489 } else {
490 upperComponentId[upperList[i]]
491 = std::vector<ttk::SimplexId>(1, upperNeighbors[i]);
492 }
493 }
494 for(auto elt : upperComponentId) {
495 upperComponents->push_back(std::vector<ttk::SimplexId>());
496 for(ttk::SimplexId i = 0; i < (ttk::SimplexId)elt.second.size(); i++) {
497 upperComponents->back().push_back(elt.second.at(i));
498 }
499 }
500
501 for(ttk::SimplexId i = 0; i < (SimplexId)lowerNeighbors.size(); i++) {
502 it = lowerComponentId.find(lowerList[i]);
503 if(it != lowerComponentId.end()) {
504 lowerComponentId[lowerList[i]].push_back(lowerNeighbors[i]);
505 } else {
506 lowerComponentId[lowerList[i]]
507 = std::vector<ttk::SimplexId>(1, lowerNeighbors[i]);
508 }
509 }
510 for(auto elt : lowerComponentId) {
511 lowerComponents->push_back(std::vector<ttk::SimplexId>());
512 for(ttk::SimplexId i = 0; i < (ttk::SimplexId)elt.second.size(); i++) {
513 lowerComponents->back().push_back(elt.second.at(i));
514 }
515 }
516
517 if(debugLevel_ >= (int)(debug::Priority::VERBOSE)) {
518 printMsg("Vertex #" + std::to_string(vertexId)
519 + ": lowerLink-#CC=" + std::to_string(lowerComponentId.size())
520 + " upperLink-#CC=" + std::to_string(upperComponentId.size()),
522 }
523
524 return 0;
525}
526
527template <class triangulationType>
529 const SimplexId &vertexId,
530 const SimplexId *const offsets,
531 const triangulationType *triangulation,
532 std::vector<std::vector<ttk::SimplexId>> *upperComponents,
533 std::vector<std::vector<ttk::SimplexId>> *lowerComponents) const {
534
535 bool isLowerOnBoundary = false, isUpperOnBoundary = false;
536 std::vector<std::vector<ttk::SimplexId>> localUpperComponents;
537 std::vector<std::vector<ttk::SimplexId>> localLowerComponents;
538 if(upperComponents == nullptr) {
539 upperComponents = &localUpperComponents;
540 }
541 if(lowerComponents == nullptr) {
542 lowerComponents = &localLowerComponents;
543 }
544 getLowerUpperComponents(vertexId, offsets, triangulation, isLowerOnBoundary,
545 isUpperOnBoundary, upperComponents, lowerComponents);
546 ttk::SimplexId const lowerComponentNumber = lowerComponents->size();
547 ttk::SimplexId const upperComponentNumber = upperComponents->size();
548
549 if(dimension_ == 1) {
550 if(lowerComponentNumber == 0 && upperComponentNumber != 0) {
551 return (char)(CriticalType::Local_minimum);
552 } else if(lowerComponentNumber != 0 && upperComponentNumber == 0) {
553 return (char)(CriticalType::Local_maximum);
554 } else if(lowerComponentNumber == 1 && upperComponentNumber == 1) {
555 return (char)(CriticalType::Regular);
556 }
557 return (char)(CriticalType::Saddle1);
558 }
559
560 if(lowerComponentNumber == 0 && upperComponentNumber == 1) {
561 return (char)(CriticalType::Local_minimum);
562 } else if(lowerComponentNumber == 1 && upperComponentNumber == 0) {
563 return (char)(CriticalType::Local_maximum);
564 } else if(lowerComponentNumber == 1 && upperComponentNumber == 1) {
565
566 if((dimension_ == 3) && (triangulation->isVertexOnBoundary(vertexId))) {
567 // special case of boundary saddles
568 if((isUpperOnBoundary) && (!isLowerOnBoundary))
569 return (char)(CriticalType::Saddle1);
570 if((!isUpperOnBoundary) && (isLowerOnBoundary))
571 return (char)(CriticalType::Saddle2);
572 }
573
574 // regular point
575 return (char)(CriticalType::Regular);
576 } else {
577 // saddles
578 if(dimension_ == 2 || dimension_ == 1) {
579 if((lowerComponentNumber == 2 && upperComponentNumber == 1)
580 || (lowerComponentNumber == 1 && upperComponentNumber == 2)
581 || (lowerComponentNumber == 2 && upperComponentNumber == 2)) {
582 // regular saddle
583 return (char)(CriticalType::Saddle1);
584 } else {
585 // monkey saddle, saddle + extremum
586 return (char)(CriticalType::Degenerate);
587 // NOTE: you may have multi-saddles on the boundary in that
588 // configuration
589 // to make this computation 100% correct, one would need to
590 // disambiguate boundary from interior vertices
591 }
592 } else if(dimension_ == 3) {
593 if(lowerComponentNumber == 2 && upperComponentNumber == 1) {
594 return (char)(CriticalType::Saddle1);
595 } else if(lowerComponentNumber == 1 && upperComponentNumber == 2) {
596 return (char)(CriticalType::Saddle2);
597 } else {
598 // monkey saddle, saddle + extremum
599 return (char)(CriticalType::Degenerate);
600 // NOTE: we may have a similar effect in 3D (TODO)
601 }
602 }
603 }
604
605 // -2: regular points
606 return (char)(CriticalType::Regular);
607}
608
609template <class triangulationType>
611 const triangulationType *ttkNotUsed(triangulation)) {
612 if(BackEnd == BACKEND::PROGRESSIVE_TOPOLOGY
613 && !std::is_same<ttk::ImplicitWithPreconditions, triangulationType>::value
614 && !std::is_same<ttk::ImplicitNoPreconditions, triangulationType>::value) {
615
617 printWrn("Explicit, Compact or Periodic");
618 printWrn("triangulation detected.");
619 printWrn("Defaulting to the generic backend.");
621
622 BackEnd = BACKEND::GENERIC;
623 }
624
625#ifdef TTK_ENABLE_MPI
626 if((ttk::hasInitializedMPI()) && (ttk::isRunningWithMPI())) {
628 printWrn("MPI run detected.");
629 printWrn("Defaulting to the generic backend.");
631
632 BackEnd = BACKEND::GENERIC;
633 }
634#endif // TTK_ENABLE_MPI
635}
#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...
virtual SimplexId getNumberOfVertices() const
virtual int getDimensionality() const
Minimalist debugging class.
Definition Debug.h:88
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
ImplicitTriangulation is a class that provides time and memory efficient traversal methods on triangu...
TTK processing package for progressive Topological Data Analysis.
TTK processing package for the computation of critical points in PL scalar fields defined on PL manif...
void checkProgressivityRequirement(const triangulationType *triangulation)
int executeProgressive(const SimplexId *const offsets, const triangulationType *triangulation)
void preconditionTriangulation(AbstractTriangulation *triangulation)
void setVertexLinkEdgeLists(const std::vector< std::vector< std::pair< SimplexId, SimplexId > > > *edgeList)
int setVertexNumber(const SimplexId &vertexNumber)
const std::vector< std::vector< std::pair< SimplexId, SimplexId > > > * vertexLinkEdgeLists_
void setOutput(std::vector< std::pair< SimplexId, char > > *criticalPoints)
int execute(const SimplexId *const offsets, const triangulationType *triangulation)
int executeLegacy(const SimplexId *const offsets, const triangulationType *triangulation)
char getCriticalType(const SimplexId &vertexId, const SimplexId *const offsets, const triangulationType *triangulation, std::vector< std::vector< ttk::SimplexId > > *upperComponents=nullptr, std::vector< std::vector< ttk::SimplexId > > *lowerComponents=nullptr) const
void setDomainDimension(const int &dimension)
int getLowerUpperComponents(const SimplexId vertexId, const SimplexId *const offsets, const triangulationType *triangulation, bool &isLowerOnBoundary, bool &isUpperOnBoundary, std::vector< std::vector< ttk::SimplexId > > *upperComponents, std::vector< std::vector< ttk::SimplexId > > *lowerComponents) const
std::vector< std::pair< SimplexId, char > > * criticalPoints_
double getElapsedTime()
Definition Timer.h:15
static UnionFind * makeUnion(UnionFind *uf0, UnionFind *uf1)
Definition UnionFind.h:40
std::string to_string(__int128)
Definition ripserpy.cpp:99
The Topology ToolKit.
COMMON_EXPORTS int MPIrank_
Definition BaseClass.cpp:9
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)