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 chunkSize = std::max(1000, (int)vertexNumber_ / (threadNumber_ * 100));
250#endif
251
252 if(triangulation) {
253
254#ifdef TTK_ENABLE_OPENMP
255#pragma omp parallel for schedule(dynamic, chunkSize) num_threads(threadNumber_)
256#endif
257 for(SimplexId i = 0; i < (SimplexId)vertexNumber_; i++) {
258#ifdef TTK_ENABLE_MPI
259 if(!isRunningWithMPI()
260 || (isRunningWithMPI()
261 && (triangulation->getVertexRank(i) == ttk::MPIrank_)))
262#endif // TTK_ENABLE_MPI
263 {
264 vertexTypes[i] = getCriticalType(i, offsets, triangulation);
265 }
266 }
267 } else if(vertexLinkEdgeLists_) {
268 // legacy implementation
269#ifdef TTK_ENABLE_OPENMP
270#pragma omp parallel for schedule(dynamic, chunkSize) num_threads(threadNumber_)
271#endif
272 for(SimplexId i = 0; i < (SimplexId)vertexNumber_; i++) {
273 // can't use MPI there since triangulation (~> rankArray) is nullptr
274 vertexTypes[i] = getCriticalType(i, offsets, (*vertexLinkEdgeLists_)[i]);
275 }
276 }
277
278 SimplexId minimumNumber = 0, maximumNumber = 0, saddleNumber = 0,
279 oneSaddleNumber = 0, twoSaddleNumber = 0, monkeySaddleNumber = 0;
280
281 // debug msg
282 if(debugLevel_ >= (int)debug::Priority::INFO) {
283 if(dimension_ == 3) {
284 for(SimplexId i = 0; i < vertexNumber_; i++) {
285 switch(vertexTypes[i]) {
286
287 case(char)(CriticalType::Local_minimum):
288 minimumNumber++;
289 break;
290
291 case(char)(CriticalType::Saddle1):
292 oneSaddleNumber++;
293 break;
294
295 case(char)(CriticalType::Saddle2):
296 twoSaddleNumber++;
297 break;
298
299 case(char)(CriticalType::Local_maximum):
300 maximumNumber++;
301 break;
302
303 case(char)(CriticalType::Degenerate):
304 monkeySaddleNumber++;
305 break;
306 }
307 }
308 } else if(dimension_ == 2) {
309 for(SimplexId i = 0; i < vertexNumber_; i++) {
310 switch(vertexTypes[i]) {
311
312 case(char)(CriticalType::Local_minimum):
313 minimumNumber++;
314 break;
315
316 case(char)(CriticalType::Saddle1):
317 saddleNumber++;
318 break;
319
320 case(char)(CriticalType::Local_maximum):
321 maximumNumber++;
322 break;
323
324 case(char)(CriticalType::Degenerate):
325 monkeySaddleNumber++;
326 break;
327 }
328 }
329 }
330
331 {
332 std::vector<std::vector<std::string>> stats;
333 stats.push_back({" #Minima", std::to_string(minimumNumber)});
334 if(dimension_ == 3) {
335 stats.push_back({" #1-saddles", std::to_string(oneSaddleNumber)});
336 stats.push_back({" #2-saddles", std::to_string(twoSaddleNumber)});
337 }
338 if(dimension_ == 2) {
339 stats.push_back({" #Saddles", std::to_string(saddleNumber)});
340 }
341 stats.push_back({" #Multi-saddles", std::to_string(monkeySaddleNumber)});
342 stats.push_back({" #Maxima", std::to_string(maximumNumber)});
343
344 printMsg(stats);
345 }
346 }
347
348 // prepare the output
349 criticalPoints_->clear();
350 criticalPoints_->reserve(vertexNumber_);
351 for(SimplexId i = 0; i < vertexNumber_; i++) {
352 if(vertexTypes[i] != (char)(CriticalType::Regular)) {
353 criticalPoints_->emplace_back(i, vertexTypes[i]);
354 }
355 }
356
357 printMsg("Processed " + std::to_string(vertexNumber_) + " vertices", 1,
358 t.getElapsedTime(), threadNumber_);
359
360 return 0;
361}
362
363template <class triangulationType>
365 const SimplexId vertexId,
366 const SimplexId *const offsets,
367 const triangulationType *triangulation,
368 bool &isLowerOnBoundary,
369 bool &isUpperOnBoundary,
370 std::vector<std::vector<ttk::SimplexId>> *upperComponents,
371 std::vector<std::vector<ttk::SimplexId>> *lowerComponents) const {
372
373 SimplexId neighborNumber = triangulation->getVertexNeighborNumber(vertexId);
374 std::vector<SimplexId> lowerNeighbors, upperNeighbors;
375
376 for(SimplexId i = 0; i < neighborNumber; i++) {
377 SimplexId neighborId = 0;
378 triangulation->getVertexNeighbor(vertexId, i, neighborId);
379
380 if(offsets[neighborId] < offsets[vertexId]) {
381 lowerNeighbors.push_back(neighborId);
382 if(dimension_ == 3) {
383 if(triangulation->isVertexOnBoundary(neighborId))
384 isLowerOnBoundary = true;
385 }
386 }
387
388 // upper link
389 if(offsets[neighborId] > offsets[vertexId]) {
390 upperNeighbors.push_back(neighborId);
391 if(dimension_ == 3) {
392 if(triangulation->isVertexOnBoundary(neighborId))
393 isUpperOnBoundary = true;
394 }
395 }
396 }
397
398 // now do the actual work
399 std::vector<UnionFind> lowerSeeds(lowerNeighbors.size());
400 std::vector<UnionFind *> lowerList(lowerNeighbors.size());
401 std::vector<UnionFind> upperSeeds(upperNeighbors.size());
402 std::vector<UnionFind *> upperList(upperNeighbors.size());
403
404 for(SimplexId i = 0; i < (SimplexId)lowerSeeds.size(); i++) {
405 lowerList[i] = &(lowerSeeds[i]);
406 }
407 for(SimplexId i = 0; i < (SimplexId)upperSeeds.size(); i++) {
408 upperList[i] = &(upperSeeds[i]);
409 }
410
411 SimplexId vertexStarSize = triangulation->getVertexStarNumber(vertexId);
412
413 for(SimplexId i = 0; i < vertexStarSize; i++) {
414 SimplexId cellId = 0;
415 triangulation->getVertexStar(vertexId, i, cellId);
416
417 SimplexId cellSize = triangulation->getCellVertexNumber(cellId);
418 for(SimplexId j = 0; j < cellSize; j++) {
419 SimplexId neighborId0 = -1;
420 triangulation->getCellVertex(cellId, j, neighborId0);
421
422 if(neighborId0 != vertexId) {
423 // we are on the link
424
425 bool lower0 = offsets[neighborId0] < offsets[vertexId];
426
427 // connect it to everybody except himself and vertexId
428 for(SimplexId k = j + 1; k < cellSize; k++) {
429
430 SimplexId neighborId1 = -1;
431 triangulation->getCellVertex(cellId, k, neighborId1);
432
433 if((neighborId1 != neighborId0) && (neighborId1 != vertexId)) {
434
435 bool lower1 = offsets[neighborId1] < offsets[vertexId];
436
437 std::vector<SimplexId> *neighbors = &lowerNeighbors;
438 std::vector<UnionFind *> *seeds = &lowerList;
439
440 if(!lower0) {
441 neighbors = &upperNeighbors;
442 seeds = &upperList;
443 }
444
445 if(lower0 == lower1) {
446 // connect their union-find sets!
447 SimplexId lowerId0 = -1, lowerId1 = -1;
448 for(SimplexId l = 0; l < (SimplexId)neighbors->size(); l++) {
449 if((*neighbors)[l] == neighborId0) {
450 lowerId0 = l;
451 }
452 if((*neighbors)[l] == neighborId1) {
453 lowerId1 = l;
454 }
455 }
456 if((lowerId0 != -1) && (lowerId1 != -1)) {
457 (*seeds)[lowerId0] = UnionFind::makeUnion(
458 (*seeds)[lowerId0], (*seeds)[lowerId1]);
459 (*seeds)[lowerId1] = (*seeds)[lowerId0];
460 }
461 }
462 }
463 }
464 }
465 }
466 }
467
468 // let's remove duplicates now
469
470 // update the UF if necessary
471 for(SimplexId i = 0; i < (SimplexId)lowerList.size(); i++)
472 lowerList[i] = lowerList[i]->find();
473 for(SimplexId i = 0; i < (SimplexId)upperList.size(); i++)
474 upperList[i] = upperList[i]->find();
475
476 std::unordered_map<UnionFind *, std::vector<ttk::SimplexId>>::iterator it;
477 std::unordered_map<UnionFind *, std::vector<ttk::SimplexId>>
478 upperComponentId{};
479 std::unordered_map<UnionFind *, std::vector<ttk::SimplexId>>
480 lowerComponentId{};
481
482 // We retrieve the lower and upper components if we want them
483 for(ttk::SimplexId i = 0; i < (SimplexId)upperNeighbors.size(); i++) {
484 it = upperComponentId.find(upperList[i]);
485 if(it != upperComponentId.end()) {
486 upperComponentId[upperList[i]].push_back(upperNeighbors[i]);
487 } else {
488 upperComponentId[upperList[i]]
489 = std::vector<ttk::SimplexId>(1, upperNeighbors[i]);
490 }
491 }
492 for(auto elt : upperComponentId) {
493 upperComponents->push_back(std::vector<ttk::SimplexId>());
494 for(ttk::SimplexId i = 0; i < (ttk::SimplexId)elt.second.size(); i++) {
495 upperComponents->back().push_back(elt.second.at(i));
496 }
497 }
498
499 for(ttk::SimplexId i = 0; i < (SimplexId)lowerNeighbors.size(); i++) {
500 it = lowerComponentId.find(lowerList[i]);
501 if(it != lowerComponentId.end()) {
502 lowerComponentId[lowerList[i]].push_back(lowerNeighbors[i]);
503 } else {
504 lowerComponentId[lowerList[i]]
505 = std::vector<ttk::SimplexId>(1, lowerNeighbors[i]);
506 }
507 }
508 for(auto elt : lowerComponentId) {
509 lowerComponents->push_back(std::vector<ttk::SimplexId>());
510 for(ttk::SimplexId i = 0; i < (ttk::SimplexId)elt.second.size(); i++) {
511 lowerComponents->back().push_back(elt.second.at(i));
512 }
513 }
514
515 if(debugLevel_ >= (int)(debug::Priority::VERBOSE)) {
516 printMsg("Vertex #" + std::to_string(vertexId)
517 + ": lowerLink-#CC=" + std::to_string(lowerComponentId.size())
518 + " upperLink-#CC=" + std::to_string(upperComponentId.size()),
520 }
521
522 return 0;
523}
524
525template <class triangulationType>
527 const SimplexId &vertexId,
528 const SimplexId *const offsets,
529 const triangulationType *triangulation,
530 std::vector<std::vector<ttk::SimplexId>> *upperComponents,
531 std::vector<std::vector<ttk::SimplexId>> *lowerComponents) const {
532
533 bool isLowerOnBoundary = false, isUpperOnBoundary = false;
534 std::vector<std::vector<ttk::SimplexId>> localUpperComponents;
535 std::vector<std::vector<ttk::SimplexId>> localLowerComponents;
536 if(upperComponents == nullptr) {
537 upperComponents = &localUpperComponents;
538 }
539 if(lowerComponents == nullptr) {
540 lowerComponents = &localLowerComponents;
541 }
542 getLowerUpperComponents(vertexId, offsets, triangulation, isLowerOnBoundary,
543 isUpperOnBoundary, upperComponents, lowerComponents);
544 ttk::SimplexId lowerComponentNumber = lowerComponents->size();
545 ttk::SimplexId upperComponentNumber = upperComponents->size();
546
547 if(dimension_ == 1) {
548 if(lowerComponentNumber == 0 && upperComponentNumber != 0) {
549 return (char)(CriticalType::Local_minimum);
550 } else if(lowerComponentNumber != 0 && upperComponentNumber == 0) {
551 return (char)(CriticalType::Local_maximum);
552 } else if(lowerComponentNumber == 1 && upperComponentNumber == 1) {
553 return (char)(CriticalType::Regular);
554 }
555 return (char)(CriticalType::Saddle1);
556 }
557
558 if(lowerComponentNumber == 0 && upperComponentNumber == 1) {
559 return (char)(CriticalType::Local_minimum);
560 } else if(lowerComponentNumber == 1 && upperComponentNumber == 0) {
561 return (char)(CriticalType::Local_maximum);
562 } else if(lowerComponentNumber == 1 && upperComponentNumber == 1) {
563
564 if((dimension_ == 3) && (triangulation->isVertexOnBoundary(vertexId))) {
565 // special case of boundary saddles
566 if((isUpperOnBoundary) && (!isLowerOnBoundary))
567 return (char)(CriticalType::Saddle1);
568 if((!isUpperOnBoundary) && (isLowerOnBoundary))
569 return (char)(CriticalType::Saddle2);
570 }
571
572 // regular point
573 return (char)(CriticalType::Regular);
574 } else {
575 // saddles
576 if(dimension_ == 2 || dimension_ == 1) {
577 if((lowerComponentNumber == 2 && upperComponentNumber == 1)
578 || (lowerComponentNumber == 1 && upperComponentNumber == 2)
579 || (lowerComponentNumber == 2 && upperComponentNumber == 2)) {
580 // regular saddle
581 return (char)(CriticalType::Saddle1);
582 } else {
583 // monkey saddle, saddle + extremum
584 return (char)(CriticalType::Degenerate);
585 // NOTE: you may have multi-saddles on the boundary in that
586 // configuration
587 // to make this computation 100% correct, one would need to
588 // disambiguate boundary from interior vertices
589 }
590 } else if(dimension_ == 3) {
591 if(lowerComponentNumber == 2 && upperComponentNumber == 1) {
592 return (char)(CriticalType::Saddle1);
593 } else if(lowerComponentNumber == 1 && upperComponentNumber == 2) {
594 return (char)(CriticalType::Saddle2);
595 } else {
596 // monkey saddle, saddle + extremum
597 return (char)(CriticalType::Degenerate);
598 // NOTE: we may have a similar effect in 3D (TODO)
599 }
600 }
601 }
602
603 // -2: regular points
604 return (char)(CriticalType::Regular);
605}
606
607template <class triangulationType>
609 const triangulationType *ttkNotUsed(triangulation)) {
610 if(BackEnd == BACKEND::PROGRESSIVE_TOPOLOGY
611 && !std::is_same<ttk::ImplicitWithPreconditions, triangulationType>::value
612 && !std::is_same<ttk::ImplicitNoPreconditions, triangulationType>::value) {
613
614 printWrn("Explicit, Compact or Periodic triangulation detected.");
615 printWrn("Defaulting to the generic backend.");
616
617 BackEnd = BACKEND::GENERIC;
618 }
619}
#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 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
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
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::GREEN+"                           "+debug::output::ENDCOLOR+debug::output::GREEN+"▒"+debug::output::ENDCOLOR+debug::output::GREEN+"▒▒▒▒▒▒▒▒▒▒▒▒▒░"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)