TTK
Loading...
Searching...
No Matches
DiscreteGradient_Template.h
Go to the documentation of this file.
1
14
15#pragma once
16
17#include <DiscreteGradient.h>
18
19using ttk::SimplexId;
20using ttk::dcg::Cell;
23
24template <typename dataType, typename triangulationType>
26 const Cell &up,
27 const Cell &down,
28 const dataType *const scalars,
29 const triangulationType &triangulation) const {
30
31 return scalars[getCellGreaterVertex(up, triangulation)]
32 - scalars[getCellGreaterVertex(down, triangulation)];
33}
34
35template <typename triangulationType>
36int DiscreteGradient::buildGradient(const triangulationType &triangulation,
37 bool bypassCache,
38 const std::vector<bool> *updateMask) {
39
40 auto &cacheHandler = *triangulation.getGradientCacheHandler();
41 const auto findGradient
42 = [this, &cacheHandler]() -> AbstractTriangulation::gradientType * {
43 if(this->inputScalarField_.first == nullptr) {
44 return {};
45 }
46 return cacheHandler.get(this->inputScalarField_);
47 };
48
49#ifdef TTK_ENABLE_OPENMP
50 if(!bypassCache && omp_in_parallel()) {
51 this->printWrn(
52 "buildGradient() called inside a parallel region, disabling cache...");
53 bypassCache = true;
54 }
55#endif // TTK_ENABLE_OPENMP
56
57 // set member variables at each buildGradient() call
58 this->dimensionality_ = triangulation.getCellVertexNumber(0) - 1;
59 this->numberOfVertices_ = triangulation.getNumberOfVertices();
60
61 this->gradient_ = bypassCache ? &this->localGradient_ : findGradient();
62 if(this->gradient_ == nullptr || bypassCache) {
63
64 if(!bypassCache) {
65 // add new cache entry
66 cacheHandler.insert(this->inputScalarField_, {});
67 this->gradient_ = cacheHandler.get(this->inputScalarField_);
68 }
69
70 // allocate gradient memory
71 this->initMemory(triangulation);
72
73 Timer tm{};
74 // compute gradient pairs
75 if(updateMask) {
76 this->processLowerStarsWithMask(
77 this->inputOffsets_, triangulation, updateMask);
78 this->printMsg("Update cached discrete gradient", 1.0,
79 tm.getElapsedTime(), this->threadNumber_);
80 } else {
81 this->processLowerStars(this->inputOffsets_, triangulation);
82 this->printMsg("Built discrete gradient", 1.0, tm.getElapsedTime(),
83 this->threadNumber_);
84 }
85 } else {
86 this->printMsg("Fetched cached discrete gradient");
87 if(updateMask) {
88 Timer tm{};
89 this->processLowerStarsWithMask(
90 this->inputOffsets_, triangulation, updateMask);
91 this->printMsg("Update cached discrete gradient", 1.0,
92 tm.getElapsedTime(), this->threadNumber_);
93 }
94 }
95
96 return 0;
97}
98
99template <typename triangulationType>
101 const std::array<std::vector<SimplexId>, 4> &criticalCellsByDim,
102 std::vector<std::array<float, 3>> &points,
103 std::vector<char> &cellDimensions,
104 std::vector<SimplexId> &cellIds,
105 std::vector<char> &isOnBoundary,
106 std::vector<SimplexId> &PLVertexIdentifiers,
107 const triangulationType &triangulation) const {
108
109 std::array<size_t, 5> partSums{};
110 for(size_t i = 0; i < criticalCellsByDim.size(); ++i) {
111 partSums[i + 1] = partSums[i] + criticalCellsByDim[i].size();
112 }
113
114 const auto nCritPoints = partSums.back();
115
116 points.resize(nCritPoints);
117 cellDimensions.resize(nCritPoints);
118 cellIds.resize(nCritPoints);
119 isOnBoundary.resize(nCritPoints);
120 PLVertexIdentifiers.resize(nCritPoints);
121
122 for(size_t i = 0; i < criticalCellsByDim.size(); ++i) {
123#ifdef TTK_ENABLE_OPENMP
124#pragma omp parallel for num_threads(threadNumber_)
125#endif // TTK_ENABLE_OPENMP
126 for(size_t j = 0; j < criticalCellsByDim[i].size(); ++j) {
127 const SimplexId cellId = criticalCellsByDim[i][j];
128 const int cellDim = i;
129 const auto o{partSums[i] + j};
130
131 triangulation.getCellIncenter(cellId, i, points[o].data());
132 cellDimensions[o] = cellDim;
133#ifdef TTK_ENABLE_MPI
134 ttk::SimplexId globalId{-1};
135 triangulation.getDistributedGlobalCellId(cellId, cellDim, globalId);
136 cellIds[o] = globalId;
137#else
138 cellIds[o] = cellId;
139#endif // TTK_ENABLE_MPI
140 const Cell cell{static_cast<int>(i), cellId};
141 isOnBoundary[o] = this->isBoundary(cell, triangulation);
142 PLVertexIdentifiers[o] = this->getCellGreaterVertex(cell, triangulation);
143 }
144 }
145
146 std::vector<std::vector<std::string>> rows(this->dimensionality_ + 1);
147 for(int i = 0; i < this->dimensionality_ + 1; ++i) {
148 rows[i]
149 = std::vector<std::string>{"#" + std::to_string(i) + "-cell(s)",
150 std::to_string(criticalCellsByDim[i].size())};
151 }
152 this->printMsg(rows);
153
154 return 0;
155}
156
157template <typename triangulationType>
159 std::vector<std::array<float, 3>> &points,
160 std::vector<char> &cellDimensions,
161 std::vector<SimplexId> &cellIds,
162 std::vector<char> &isOnBoundary,
163 std::vector<SimplexId> &PLVertexIdentifiers,
164 const triangulationType &triangulation) const {
165
166 std::array<std::vector<SimplexId>, 4> criticalCellsByDim;
167 getCriticalPoints(criticalCellsByDim, triangulation);
168 setCriticalPoints(criticalCellsByDim, points, cellDimensions, cellIds,
169 isOnBoundary, PLVertexIdentifiers, triangulation);
170
171 return 0;
172}
173
174template <typename triangulationType>
176 std::array<std::vector<SimplexId>, 4> &criticalCellsByDim,
177 const triangulationType &triangulation) const {
178
179 const auto dims{this->getNumberOfDimensions()};
180 for(int i = 0; i < dims; ++i) {
181
182 // map: store critical cell per dimension per thread
183 std::vector<std::vector<SimplexId>> critCellsPerThread(this->threadNumber_);
184 const auto numberOfCells{this->getNumberOfCells(i, triangulation)};
185
186 // use static scheduling to ensure that critical cells
187 // are sorted by id
188
189#ifdef TTK_ENABLE_OPENMP
190#pragma omp parallel for num_threads(this->threadNumber_) schedule(static)
191#endif // TTK_ENABLE_OPENMP
192 for(SimplexId j = 0; j < numberOfCells; ++j) {
193#ifdef TTK_ENABLE_OPENMP
194 const auto tid = omp_get_thread_num();
195#else
196 const auto tid = 0;
197#endif // TTK_ENABLE_OPENMP
198 if(this->isCellCritical(i, j)) {
199 critCellsPerThread[tid].emplace_back(j);
200 }
201 }
202
203 // reduce: aggregate critical cells per thread
204 criticalCellsByDim[i] = std::move(critCellsPerThread[0]);
205 for(size_t j = 1; j < critCellsPerThread.size(); ++j) {
206 const auto &vec{critCellsPerThread[j]};
207 criticalCellsByDim[i].insert(
208 criticalCellsByDim[i].end(), vec.begin(), vec.end());
209 }
210 }
211
212 return 0;
213}
214
215template <typename triangulationType>
217 const int dimension, const triangulationType &triangulation) const {
218
219 if(dimension > this->dimensionality_ || dimension < 0) {
220 return -1;
221 }
222
223 switch(dimension) {
224 case 0:
225 return triangulation.getNumberOfVertices();
226 break;
227
228 case 1:
229 return triangulation.getNumberOfEdges();
230 break;
231
232 case 2:
233 return triangulation.getNumberOfTriangles();
234 break;
235
236 case 3:
237 return triangulation.getNumberOfCells();
238 break;
239 }
240
241 return -1;
242}
243
244template <typename triangulationType>
245inline void DiscreteGradient::lowerStarWithMask(
246 lowerStarType &ls,
247 const SimplexId a,
248 const SimplexId *const offsets,
249 const triangulationType &triangulation,
250 const std::vector<bool> *updateMask) const {
251
252 // make sure that ls is cleared
253 for(auto &vec : ls) {
254 vec.clear();
255 }
256
257 // a belongs to its lower star
258 ls[0].emplace_back(CellExt{0, a});
259
260 if(updateMask != nullptr) {
261 for(auto &vec : ls[0]) {
262 int vertexId = vec.id_;
263 int edgeId = (*gradient_)[0][vertexId];
264
265 (*gradient_)[0][vertexId] = -1;
266 if(edgeId != -1) {
267 (*gradient_)[1][edgeId] = -1;
268 }
269 }
270 }
271
272 // store lower edges
273 const auto nedges = triangulation.getVertexEdgeNumber(a);
274 ls[1].reserve(nedges);
275 for(SimplexId i = 0; i < nedges; i++) {
276 SimplexId edgeId;
277 triangulation.getVertexEdge(a, i, edgeId);
278 SimplexId vertexId;
279 triangulation.getEdgeVertex(edgeId, 0, vertexId);
280 if(vertexId == a) {
281 triangulation.getEdgeVertex(edgeId, 1, vertexId);
282 }
283 if(offsets[vertexId] < offsets[a]) {
284 ls[1].emplace_back(CellExt{1, edgeId, {offsets[vertexId], -1, -1}, {}});
285 }
286 }
287
288 if(updateMask != nullptr) {
289 for(auto &vec : ls[1]) {
290 int edgeId = vec.id_;
291 int vertexId = (*gradient_)[1][edgeId];
292 int triangleId = (*gradient_)[2][edgeId];
293
294 (*gradient_)[1][edgeId] = -1;
295 if(vertexId != -1) {
296 (*gradient_)[0][vertexId] = -1;
297 }
298
299 (*gradient_)[2][edgeId] = -1;
300 if(triangleId != -1) {
301 (*gradient_)[3][triangleId] = -1;
302 }
303 }
304 }
305
306 if(ls[1].size() < 2) {
307 // at least two edges in the lower star for one triangle
308 return;
309 }
310
311 const auto processTriangle
312 = [&](const SimplexId triangleId, const SimplexId v0, const SimplexId v1,
313 const SimplexId v2) {
314 std::array<SimplexId, 3> lowVerts{-1, -1, -1};
315 if(v0 == a) {
316 lowVerts[0] = offsets[v1];
317 lowVerts[1] = offsets[v2];
318 } else if(v1 == a) {
319 lowVerts[0] = offsets[v0];
320 lowVerts[1] = offsets[v2];
321 } else if(v2 == a) {
322 lowVerts[0] = offsets[v0];
323 lowVerts[1] = offsets[v1];
324 }
325 // higher order vertex first
326 if(lowVerts[0] < lowVerts[1]) {
327 std::swap(lowVerts[0], lowVerts[1]);
328 }
329 if(offsets[a] > lowVerts[0]) { // triangle in lowerStar
330 uint8_t j{}, k{};
331 // store edges indices of current triangle
332 std::array<uint8_t, 3> faces{};
333 for(const auto &e : ls[1]) {
334 if(e.lowVerts_[0] == lowVerts[0] || e.lowVerts_[0] == lowVerts[1]) {
335 faces[k++] = j;
336 }
337 j++;
338 }
339 ls[2].emplace_back(CellExt{2, triangleId, lowVerts, faces});
340 }
341 };
342
343 if(dimensionality_ == 2) {
344 // store lower triangles
345
346 // use optimised triangulation methods:
347 // getVertexStar instead of getVertexTriangle
348 // getCellVertex instead of getTriangleVertex
349 const auto ncells = triangulation.getVertexStarNumber(a);
350 ls[2].reserve(ncells);
351 for(SimplexId i = 0; i < ncells; ++i) {
352 SimplexId cellId;
353 triangulation.getVertexStar(a, i, cellId);
354 SimplexId v0{}, v1{}, v2{};
355 triangulation.getCellVertex(cellId, 0, v0);
356 triangulation.getCellVertex(cellId, 1, v1);
357 triangulation.getCellVertex(cellId, 2, v2);
358 processTriangle(cellId, v0, v1, v2);
359 }
360
361 if(updateMask != nullptr) {
362 for(auto &vec : ls[2]) {
363 int triangleId = vec.id_;
364 int edgeId = (*gradient_)[3][triangleId];
365
366 (*gradient_)[3][triangleId] = -1;
367 if(edgeId != -1) {
368 (*gradient_)[2][edgeId] = -1;
369 }
370 }
371 }
372
373 } else if(dimensionality_ == 3) {
374 // store lower triangles
375 const auto ntri = triangulation.getVertexTriangleNumber(a);
376 ls[2].reserve(ntri);
377 for(SimplexId i = 0; i < ntri; i++) {
378 SimplexId triangleId;
379 triangulation.getVertexTriangle(a, i, triangleId);
380 SimplexId v0{}, v1{}, v2{};
381 triangulation.getTriangleVertex(triangleId, 0, v0);
382 triangulation.getTriangleVertex(triangleId, 1, v1);
383 triangulation.getTriangleVertex(triangleId, 2, v2);
384 processTriangle(triangleId, v0, v1, v2);
385 }
386
387 if(updateMask != nullptr) {
388 for(auto &vec : ls[2]) {
389 int triangleId = vec.id_;
390 int edgeId = (*gradient_)[3][triangleId];
391 int tetraId = (*gradient_)[4][triangleId];
392
393 (*gradient_)[3][triangleId] = -1;
394 if(edgeId != -1) {
395 (*gradient_)[2][edgeId] = -1;
396 }
397
398 (*gradient_)[4][triangleId] = -1;
399 if(tetraId != -1) {
400 (*gradient_)[5][tetraId] = -1;
401 }
402 }
403 }
404
405 // at least three triangles in the lower star for one tetra
406 if(ls[2].size() >= 3) {
407 // store lower tetra
408 const auto ncells = triangulation.getVertexStarNumber(a);
409 ls[3].reserve(ncells);
410 for(SimplexId i = 0; i < ncells; ++i) {
411 SimplexId cellId;
412 triangulation.getVertexStar(a, i, cellId);
413 std::array<SimplexId, 3> lowVerts{-1, -1, -1};
414 SimplexId v0{}, v1{}, v2{}, v3{};
415 triangulation.getCellVertex(cellId, 0, v0);
416 triangulation.getCellVertex(cellId, 1, v1);
417 triangulation.getCellVertex(cellId, 2, v2);
418 triangulation.getCellVertex(cellId, 3, v3);
419 if(v0 == a) {
420 lowVerts[0] = offsets[v1];
421 lowVerts[1] = offsets[v2];
422 lowVerts[2] = offsets[v3];
423 } else if(v1 == a) {
424 lowVerts[0] = offsets[v0];
425 lowVerts[1] = offsets[v2];
426 lowVerts[2] = offsets[v3];
427 } else if(v2 == a) {
428 lowVerts[0] = offsets[v0];
429 lowVerts[1] = offsets[v1];
430 lowVerts[2] = offsets[v3];
431 } else if(v3 == a) {
432 lowVerts[0] = offsets[v0];
433 lowVerts[1] = offsets[v1];
434 lowVerts[2] = offsets[v2];
435 }
436 if(offsets[a] > *std::max_element(
437 lowVerts.begin(), lowVerts.end())) { // tetra in lowerStar
438
439 // higher order vertex first
440 std::sort(lowVerts.rbegin(), lowVerts.rend());
441
442 uint8_t j{}, k{};
443 // store triangles indices of current tetra
444 std::array<uint8_t, 3> faces{};
445 for(const auto &t : ls[2]) {
446 // lowVerts & t.lowVerts are ordered, no need to check if
447 // t.lowVerts[0] == lowVerts[2] or t.lowVerts[1] == lowVerts[0]
448 if((t.lowVerts_[0] == lowVerts[0]
449 && (t.lowVerts_[1] == lowVerts[1]
450 || t.lowVerts_[1] == lowVerts[2]))
451 || (t.lowVerts_[0] == lowVerts[1]
452 && t.lowVerts_[1] == lowVerts[2])) {
453 faces[k++] = j;
454 }
455 j++;
456 }
457
458 ls[3].emplace_back(CellExt{3, cellId, lowVerts, faces});
459 }
460 }
461
462 if(updateMask != nullptr) {
463 for(auto &vec : ls[3]) {
464 int tetraId = vec.id_;
465 int triangleId = (*gradient_)[5][tetraId];
466 (*gradient_)[5][tetraId] = -1;
467 if(triangleId != -1) {
468 (*gradient_)[4][triangleId] = -1;
469 }
470 }
471 }
472 }
473 }
474}
475
476template <typename triangulationType>
477inline void
478 DiscreteGradient::lowerStar(lowerStarType &ls,
479 const SimplexId a,
480 const SimplexId *const offsets,
481 const triangulationType &triangulation) const {
482
483 // WARNING
484 // If you modify this function, please make sure to also report your edit to
485 // the other implementation of this function, `lowerStarWithMask`.
486
487 // make sure that ls is cleared
488 for(auto &vec : ls) {
489 vec.clear();
490 }
491
492 // a belongs to its lower star
493 CellExt const localCellExt{0, a};
494 ls[0].emplace_back(localCellExt);
495
496 // store lower edges
497 const auto nedges = triangulation.getVertexEdgeNumber(a);
498 ls[1].reserve(nedges);
499 for(SimplexId i = 0; i < nedges; i++) {
500 SimplexId edgeId;
501 triangulation.getVertexEdge(a, i, edgeId);
502 SimplexId vertexId;
503 triangulation.getEdgeVertex(edgeId, 0, vertexId);
504 if(vertexId == a) {
505 triangulation.getEdgeVertex(edgeId, 1, vertexId);
506 }
507 if(offsets[vertexId] < offsets[a]) {
508 ls[1].emplace_back(CellExt{1, edgeId, {offsets[vertexId], -1, -1}, {}});
509 }
510 }
511
512 if(ls[1].size() < 2) {
513 // at least two edges in the lower star for one triangle
514 return;
515 }
516
517 const auto processTriangle
518 = [&](const SimplexId triangleId, const SimplexId v0, const SimplexId v1,
519 const SimplexId v2) {
520 std::array<SimplexId, 3> lowVerts{-1, -1, -1};
521 if(v0 == a) {
522 lowVerts[0] = offsets[v1];
523 lowVerts[1] = offsets[v2];
524 } else if(v1 == a) {
525 lowVerts[0] = offsets[v0];
526 lowVerts[1] = offsets[v2];
527 } else if(v2 == a) {
528 lowVerts[0] = offsets[v0];
529 lowVerts[1] = offsets[v1];
530 }
531 // higher order vertex first
532 if(lowVerts[0] < lowVerts[1]) {
533 std::swap(lowVerts[0], lowVerts[1]);
534 }
535 if(offsets[a] > lowVerts[0]) { // triangle in lowerStar
536 uint8_t j{}, k{};
537 // store edges indices of current triangle
538 std::array<uint8_t, 3> faces{};
539 for(const auto &e : ls[1]) {
540 if(e.lowVerts_[0] == lowVerts[0] || e.lowVerts_[0] == lowVerts[1]) {
541 faces[k++] = j;
542 }
543 j++;
544 }
545 CellExt const localCellExt2{2, triangleId, lowVerts, faces};
546 ls[2].emplace_back(localCellExt2);
547 }
548 };
549
550 if(dimensionality_ == 2) {
551 // store lower triangles
552
553 // use optimised triangulation methods:
554 // getVertexStar instead of getVertexTriangle
555 // getCellVertex instead of getTriangleVertex
556 const auto ncells = triangulation.getVertexStarNumber(a);
557 ls[2].reserve(ncells);
558 for(SimplexId i = 0; i < ncells; ++i) {
559 SimplexId cellId;
560 triangulation.getVertexStar(a, i, cellId);
561 SimplexId v0{}, v1{}, v2{};
562 triangulation.getCellVertex(cellId, 0, v0);
563 triangulation.getCellVertex(cellId, 1, v1);
564 triangulation.getCellVertex(cellId, 2, v2);
565 processTriangle(cellId, v0, v1, v2);
566 }
567 } else if(dimensionality_ == 3) {
568 // store lower triangles
569 const auto ntri = triangulation.getVertexTriangleNumber(a);
570 ls[2].reserve(ntri);
571 for(SimplexId i = 0; i < ntri; i++) {
572 SimplexId triangleId;
573 triangulation.getVertexTriangle(a, i, triangleId);
574 SimplexId v0{}, v1{}, v2{};
575 triangulation.getTriangleVertex(triangleId, 0, v0);
576 triangulation.getTriangleVertex(triangleId, 1, v1);
577 triangulation.getTriangleVertex(triangleId, 2, v2);
578 processTriangle(triangleId, v0, v1, v2);
579 }
580
581 // at least three triangles in the lower star for one tetra
582 if(ls[2].size() >= 3) {
583 // store lower tetra
584 const auto ncells = triangulation.getVertexStarNumber(a);
585 ls[3].reserve(ncells);
586 for(SimplexId i = 0; i < ncells; ++i) {
587 SimplexId cellId;
588 triangulation.getVertexStar(a, i, cellId);
589 std::array<SimplexId, 3> lowVerts{-1, -1, -1};
590 SimplexId v0{}, v1{}, v2{}, v3{};
591 triangulation.getCellVertex(cellId, 0, v0);
592 triangulation.getCellVertex(cellId, 1, v1);
593 triangulation.getCellVertex(cellId, 2, v2);
594 triangulation.getCellVertex(cellId, 3, v3);
595 if(v0 == a) {
596 lowVerts[0] = offsets[v1];
597 lowVerts[1] = offsets[v2];
598 lowVerts[2] = offsets[v3];
599 } else if(v1 == a) {
600 lowVerts[0] = offsets[v0];
601 lowVerts[1] = offsets[v2];
602 lowVerts[2] = offsets[v3];
603 } else if(v2 == a) {
604 lowVerts[0] = offsets[v0];
605 lowVerts[1] = offsets[v1];
606 lowVerts[2] = offsets[v3];
607 } else if(v3 == a) {
608 lowVerts[0] = offsets[v0];
609 lowVerts[1] = offsets[v1];
610 lowVerts[2] = offsets[v2];
611 }
612 if(offsets[a] > *std::max_element(
613 lowVerts.begin(), lowVerts.end())) { // tetra in lowerStar
614
615 // higher order vertex first
616 std::sort(lowVerts.rbegin(), lowVerts.rend());
617
618 uint8_t j{}, k{};
619 // store triangles indices of current tetra
620 std::array<uint8_t, 3> faces{};
621 for(const auto &t : ls[2]) {
622 // lowVerts & t.lowVerts are ordered, no need to check if
623 // t.lowVerts[0] == lowVerts[2] or t.lowVerts[1] == lowVerts[0]
624 if((t.lowVerts_[0] == lowVerts[0]
625 && (t.lowVerts_[1] == lowVerts[1]
626 || t.lowVerts_[1] == lowVerts[2]))
627 || (t.lowVerts_[0] == lowVerts[1]
628 && t.lowVerts_[1] == lowVerts[2])) {
629 faces[k++] = j;
630 }
631 j++;
632 }
633
634 CellExt const localCellExt3{3, cellId, lowVerts, faces};
635 ls[3].emplace_back(localCellExt3);
636 }
637 }
638 }
639 }
640}
641
642template <typename triangulationType>
643inline void DiscreteGradient::pairCells(
644 CellExt &alpha, CellExt &beta, const triangulationType &triangulation) {
645#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
646 char localBId{0}, localAId{0};
647 SimplexId a{}, b{};
648
649 if(beta.dim_ == 1) {
650
651 for(SimplexId i = 0; i < 2; ++i) {
652 triangulation.getEdgeVertex(beta.id_, i, a);
653 if(a == alpha.id_) {
654 localAId = i;
655 break;
656 }
657 }
658 const auto nedges = triangulation.getVertexEdgeNumber(alpha.id_);
659 for(SimplexId i = 0; i < nedges; ++i) {
660 triangulation.getVertexEdge(alpha.id_, i, b);
661 if(b == beta.id_) {
662 localBId = i;
663 break;
664 }
665 }
666 } else if(beta.dim_ == 2) {
667 for(SimplexId i = 0; i < 3; ++i) {
668 triangulation.getTriangleEdge(beta.id_, i, a);
669 if(a == alpha.id_) {
670 localAId = i;
671 break;
672 }
673 }
674 const auto ntri = triangulation.getEdgeTriangleNumber(alpha.id_);
675 for(SimplexId i = 0; i < ntri; ++i) {
676 triangulation.getEdgeTriangle(alpha.id_, i, b);
677 if(b == beta.id_) {
678 localBId = i;
679 break;
680 }
681 }
682 } else {
683 for(SimplexId i = 0; i < 4; ++i) {
684 triangulation.getCellTriangle(beta.id_, i, a);
685 if(a == alpha.id_) {
686 localAId = i;
687 break;
688 }
689 }
690 const auto ntetra = triangulation.getTriangleStarNumber(alpha.id_);
691 for(SimplexId i = 0; i < ntetra; ++i) {
692 triangulation.getTriangleStar(alpha.id_, i, b);
693 if(b == beta.id_) {
694 localBId = i;
695 break;
696 }
697 }
698 }
699 (*gradient_)[2 * alpha.dim_][alpha.id_] = localBId;
700 (*gradient_)[2 * alpha.dim_ + 1][beta.id_] = localAId;
701#else
702 TTK_FORCE_USE(triangulation);
703 (*gradient_)[2 * alpha.dim_][alpha.id_] = beta.id_;
704 (*gradient_)[2 * alpha.dim_ + 1][beta.id_] = alpha.id_;
705#endif // TTK_ENABLE_DCG_OPTIMIZE_MEMORY
706 alpha.paired_ = true;
707 beta.paired_ = true;
708}
709
710template <typename triangulationType>
711int DiscreteGradient::processLowerStars(
712 const SimplexId *const offsets, const triangulationType &triangulation) {
713
714 // WARNING
715 // If you modify this function, please make sure to also report your edit to
716 // the other implementation of this function, `processLowerStarsWithMask`.
717
718 /* Compute gradient */
719
720 auto nverts = triangulation.getNumberOfVertices();
721
722 // Comparison function for Cells inside priority queues
723 const auto orderCells = [&](const CellExt &a, const CellExt &b) -> bool {
724 return a.lowVerts_ > b.lowVerts_;
725 };
726
727 // Type alias for priority queues
728 using pqType
729 = std::priority_queue<std::reference_wrapper<CellExt>,
730 std::vector<std::reference_wrapper<CellExt>>,
731 decltype(orderCells)>;
732
733 // To reduce allocations, priority queues and lowerStar objects are
734 // cleaned & reused between iterations.
735
736 // Priority queues are pushed at the beginning and popped at the
737 // end. To pop the minimum, elements should be sorted in a
738 // decreasing order.
739 pqType pqZero{orderCells}, pqOne{orderCells};
740
741 // store lower star structure
742 lowerStarType Lx;
743
744#ifdef TTK_ENABLE_OPENMP
745#pragma omp parallel for num_threads(threadNumber_) \
746 firstprivate(Lx, pqZero, pqOne)
747#endif // TTK_ENABLE_OPENMP
748 for(SimplexId x = 0; x < nverts; x++) {
749
750 // clear priority queues (they should be empty at the end of the
751 // previous iteration)
752 while(!pqZero.empty()) {
753 pqZero.pop();
754 }
755 while(!pqOne.empty()) {
756 pqOne.pop();
757 }
758
759 // Insert into pqOne cofacets of cell c_alpha such as numUnpairedFaces == 1
760 const auto insertCofacets = [&](const CellExt &ca, lowerStarType &ls) {
761 if(ca.dim_ == 1) {
762 for(auto &beta : ls[2]) {
763 if(ls[1][beta.faces_[0]].id_ == ca.id_
764 || ls[1][beta.faces_[1]].id_ == ca.id_) {
765 // edge ca belongs to triangle beta
766 if(numUnpairedFacesTriangle(beta, ls).first == 1) {
767 pqOne.push(beta);
768 }
769 }
770 }
771
772 } else if(ca.dim_ == 2) {
773 for(auto &beta : ls[3]) {
774 if(ls[2][beta.faces_[0]].id_ == ca.id_
775 || ls[2][beta.faces_[1]].id_ == ca.id_
776 || ls[2][beta.faces_[2]].id_ == ca.id_) {
777 // triangle ca belongs to tetra beta
778 if(numUnpairedFacesTetra(beta, ls).first == 1) {
779 pqOne.push(beta);
780 }
781 }
782 }
783 }
784 };
785
786 lowerStar(Lx, x, offsets, triangulation);
787 // In case the vertex is a ghost, the gradient of the
788 // simplices of its star is set to GHOST_GRADIENT
789#ifdef TTK_ENABLE_MPI
790 if(ttk::isRunningWithMPI()
791 && triangulation.getVertexRank(x) != ttk::MPIrank_) {
792 int sizeDim = Lx.size();
793 for(int i = 0; i < sizeDim; i++) {
794 int nCells = Lx[i].size();
795 for(int j = 0; j < nCells; j++) {
796 setCellToGhost(Lx[i][j].dim_, Lx[i][j].id_);
797 }
798 }
799 } else
800#endif // TTK_ENABLE_MPI
801
802 {
803 // Lx[1] empty => x is a local minimum
804 if(!Lx[1].empty()) {
805 // get delta: 1-cell (edge) with minimal G value (steeper gradient)
806 size_t minId = 0;
807 for(size_t i = 1; i < Lx[1].size(); ++i) {
808 const auto &a = Lx[1][minId].lowVerts_[0];
809 const auto &b = Lx[1][i].lowVerts_[0];
810 if(a > b) {
811 // edge[i] < edge[0]
812 minId = i;
813 }
814 }
815
816 auto &c_delta = Lx[1][minId];
817
818 // store x (0-cell) -> delta (1-cell) V-path
819 pairCells(Lx[0][0], c_delta, triangulation);
820
821 // push every 1-cell in Lx that is not delta into pqZero
822 for(auto &alpha : Lx[1]) {
823 if(alpha.id_ != c_delta.id_) {
824 pqZero.push(alpha);
825 }
826 }
827
828 // push into pqOne every coface of delta in Lx (2-cells only,
829 // 3-cells have not any facet paired yet) such that
830 // numUnpairedFaces == 1
831 insertCofacets(c_delta, Lx);
832
833 while(!pqOne.empty() || !pqZero.empty()) {
834 while(!pqOne.empty()) {
835 auto &c_alpha = pqOne.top().get();
836 pqOne.pop();
837 auto unpairedFaces = numUnpairedFaces(c_alpha, Lx);
838 if(unpairedFaces.first == 0) {
839 pqZero.push(c_alpha);
840 } else {
841 auto &c_pair_alpha = Lx[c_alpha.dim_ - 1][unpairedFaces.second];
842
843 // store (pair_alpha) -> (alpha) V-path
844 pairCells(c_pair_alpha, c_alpha, triangulation);
845
846 // add cofaces of c_alpha and c_pair_alpha to pqOne
847 insertCofacets(c_alpha, Lx);
848 insertCofacets(c_pair_alpha, Lx);
849 }
850 }
851
852 // skip pair_alpha from pqZero:
853 // cells in pqZero are not critical if already paired
854 while(!pqZero.empty() && pqZero.top().get().paired_) {
855 pqZero.pop();
856 }
857
858 if(!pqZero.empty()) {
859 auto &c_gamma = pqZero.top().get();
860 pqZero.pop();
861
862 // gamma is a critical cell
863 // mark gamma as paired
864 c_gamma.paired_ = true;
865
866 // add cofacets of c_gamma to pqOne
867 insertCofacets(c_gamma, Lx);
868 }
869 }
870 }
871 }
872 }
873
874 return 0;
875}
876
877template <typename triangulationType>
878int DiscreteGradient::processLowerStarsWithMask(
879 const SimplexId *const offsets,
880 const triangulationType &triangulation,
881 const std::vector<bool> *updateMask) {
882
883 auto nverts = triangulation.getNumberOfVertices();
884
885 /* Compute gradient */
886
887 // Comparison function for Cells inside priority queues
888 const auto orderCells = [&](const CellExt &a, const CellExt &b) -> bool {
889 return a.lowVerts_ > b.lowVerts_;
890 };
891
892 // Type alias for priority queues
893 using pqType
894 = std::priority_queue<std::reference_wrapper<CellExt>,
895 std::vector<std::reference_wrapper<CellExt>>,
896 decltype(orderCells)>;
897
898 // To reduce allocations, priority queues and lowerStar objects are
899 // cleaned & reused between iterations.
900
901 // Priority queues are pushed at the beginning and popped at the
902 // end. To pop the minimum, elements should be sorted in a
903 // decreasing order.
904 pqType pqZero{orderCells}, pqOne{orderCells};
905
906 // store lower star structure
907 lowerStarType Lx;
908
909#ifdef TTK_ENABLE_OPENMP
910#pragma omp parallel for num_threads(threadNumber_) \
911 firstprivate(Lx, pqZero, pqOne)
912#endif // TTK_ENABLE_OPENMP
913 for(SimplexId x = 0; x < nverts; x++) {
914
915 if((updateMask != nullptr) && ((*updateMask)[x] == false)) {
916 continue;
917 }
918
919 // clear priority queues (they should be empty at the end of the
920 // previous iteration)
921 while(!pqZero.empty()) {
922 pqZero.pop();
923 }
924 while(!pqOne.empty()) {
925 pqOne.pop();
926 }
927
928 // Insert into pqOne cofacets of cell c_alpha such as numUnpairedFaces == 1
929 const auto insertCofacets = [&](const CellExt &ca, lowerStarType &ls) {
930 if(ca.dim_ == 1) {
931 for(auto &beta : ls[2]) {
932 if(ls[1][beta.faces_[0]].id_ == ca.id_
933 || ls[1][beta.faces_[1]].id_ == ca.id_) {
934 // edge ca belongs to triangle beta
935 if(numUnpairedFacesTriangle(beta, ls).first == 1) {
936 pqOne.push(beta);
937 }
938 }
939 }
940
941 } else if(ca.dim_ == 2) {
942 for(auto &beta : ls[3]) {
943 if(ls[2][beta.faces_[0]].id_ == ca.id_
944 || ls[2][beta.faces_[1]].id_ == ca.id_
945 || ls[2][beta.faces_[2]].id_ == ca.id_) {
946 // triangle ca belongs to tetra beta
947 if(numUnpairedFacesTetra(beta, ls).first == 1) {
948 pqOne.push(beta);
949 }
950 }
951 }
952 }
953 };
954
955 lowerStarWithMask(Lx, x, offsets, triangulation, updateMask);
956 // In case the vertex is a ghost, the gradient of the
957 // simplices of its star is set to GHOST_GRADIENT
958#ifdef TTK_ENABLE_MPI
959 if(ttk::isRunningWithMPI()
960 && triangulation.getVertexRank(x) != ttk::MPIrank_) {
961 int sizeDim = Lx.size();
962 for(int i = 0; i < sizeDim; i++) {
963 int nCells = Lx[i].size();
964 for(int j = 0; j < nCells; j++) {
965 setCellToGhost(Lx[i][j].dim_, Lx[i][j].id_);
966 }
967 }
968 } else
969#endif // TTK_ENABLE_MPI
970
971 {
972 // Lx[1] empty => x is a local minimum
973 if(!Lx[1].empty()) {
974 // get delta: 1-cell (edge) with minimal G value (steeper gradient)
975 size_t minId = 0;
976 for(size_t i = 1; i < Lx[1].size(); ++i) {
977 const auto &a = Lx[1][minId].lowVerts_[0];
978 const auto &b = Lx[1][i].lowVerts_[0];
979 if(a > b) {
980 // edge[i] < edge[0]
981 minId = i;
982 }
983 }
984
985 auto &c_delta = Lx[1][minId];
986
987 // store x (0-cell) -> delta (1-cell) V-path
988 pairCells(Lx[0][0], c_delta, triangulation);
989
990 // push every 1-cell in Lx that is not delta into pqZero
991 for(auto &alpha : Lx[1]) {
992 if(alpha.id_ != c_delta.id_) {
993 pqZero.push(alpha);
994 }
995 }
996
997 // push into pqOne every coface of delta in Lx (2-cells only,
998 // 3-cells have not any facet paired yet) such that
999 // numUnpairedFaces == 1
1000 insertCofacets(c_delta, Lx);
1001
1002 while(!pqOne.empty() || !pqZero.empty()) {
1003 while(!pqOne.empty()) {
1004 auto &c_alpha = pqOne.top().get();
1005 pqOne.pop();
1006 auto unpairedFaces = numUnpairedFaces(c_alpha, Lx);
1007 if(unpairedFaces.first == 0) {
1008 pqZero.push(c_alpha);
1009 } else {
1010 auto &c_pair_alpha = Lx[c_alpha.dim_ - 1][unpairedFaces.second];
1011
1012 // store (pair_alpha) -> (alpha) V-path
1013 pairCells(c_pair_alpha, c_alpha, triangulation);
1014
1015 // add cofaces of c_alpha and c_pair_alpha to pqOne
1016 insertCofacets(c_alpha, Lx);
1017 insertCofacets(c_pair_alpha, Lx);
1018 }
1019 }
1020
1021 // skip pair_alpha from pqZero:
1022 // cells in pqZero are not critical if already paired
1023 while(!pqZero.empty() && pqZero.top().get().paired_) {
1024 pqZero.pop();
1025 }
1026
1027 if(!pqZero.empty()) {
1028 auto &c_gamma = pqZero.top().get();
1029 pqZero.pop();
1030
1031 // gamma is a critical cell
1032 // mark gamma as paired
1033 c_gamma.paired_ = true;
1034
1035 // add cofacets of c_gamma to pqOne
1036 insertCofacets(c_gamma, Lx);
1037 }
1038 }
1039 }
1040 }
1041 }
1042
1043 return 0;
1044}
1045
1046template <typename triangulationType>
1048 const Cell &cell, const triangulationType &triangulation) const {
1049
1050 if(cell.dim_ > this->dimensionality_ || cell.dim_ < 0) {
1051 return false;
1052 }
1053
1054 const auto vert{this->getCellGreaterVertex(cell, triangulation)};
1055 return triangulation.isVertexOnBoundary(vert);
1056}
1057
1058template <typename triangulationType>
1061 const triangulationType &triangulation,
1062 bool isReverse) const {
1063
1064 // ensure that getPairedCell(Cell, boolean) calls are rejected
1065 static_assert(
1066 std::is_base_of<AbstractTriangulation, triangulationType>(),
1067 "triangulationType should be an AbstractTriangulation derivative");
1068
1069#ifndef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1070 TTK_FORCE_USE(triangulation);
1071#endif // !TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1072
1073 if((cell.dim_ > this->dimensionality_ - 1 && !isReverse)
1074 || (cell.dim_ > this->dimensionality_ && isReverse) || cell.dim_ < 0) {
1075 return -1;
1076 }
1077
1078 SimplexId id{-1};
1079
1080 if(cell.dim_ == 0) {
1081 if(!isReverse) {
1082#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1083 const auto locId{(*gradient_)[0][cell.id_]};
1084 if(locId != -1) {
1085 triangulation.getVertexEdge(cell.id_, locId, id);
1086 }
1087#else
1088 id = (*gradient_)[0][cell.id_];
1089#endif
1090 }
1091 }
1092
1093 else if(cell.dim_ == 1) {
1094 if(isReverse) {
1095#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1096 const auto locId{(*gradient_)[1][cell.id_]};
1097 if(locId != -1) {
1098 triangulation.getEdgeVertex(cell.id_, locId, id);
1099 }
1100#else
1101 id = (*gradient_)[1][cell.id_];
1102#endif
1103 } else {
1104#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1105 const auto locId{(*gradient_)[2][cell.id_]};
1106 if(locId != -1) {
1107 triangulation.getEdgeTriangle(cell.id_, locId, id);
1108 }
1109#else
1110 id = (*gradient_)[2][cell.id_];
1111#endif
1112 }
1113 }
1114
1115 else if(cell.dim_ == 2) {
1116 if(isReverse) {
1117#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1118 const auto locId{(*gradient_)[3][cell.id_]};
1119 if(locId != -1) {
1120 triangulation.getTriangleEdge(cell.id_, locId, id);
1121 }
1122#else
1123 id = (*gradient_)[3][cell.id_];
1124#endif
1125 } else {
1126#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1127 const auto locId{(*gradient_)[4][cell.id_]};
1128 if(locId != -1) {
1129 triangulation.getTriangleStar(cell.id_, locId, id);
1130 }
1131#else
1132 id = (*gradient_)[4][cell.id_];
1133#endif
1134 }
1135 }
1136
1137 else if(cell.dim_ == 3) {
1138 if(isReverse) {
1139#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1140 const auto locId{(*gradient_)[5][cell.id_]};
1141 if(locId != -1) {
1142 triangulation.getCellTriangle(cell.id_, locId, id);
1143 }
1144#else
1145 id = (*gradient_)[5][cell.id_];
1146#endif
1147 }
1148 }
1149
1150 return id;
1151}
1152
1153template <typename triangulationType>
1155 const Cell &cell,
1156 std::vector<Cell> &vpath,
1157 const triangulationType &triangulation) const {
1158
1159 if(cell.dim_ == 0) {
1160 // assume that cellId is a vertex
1161 SimplexId currentId = cell.id_;
1162 SimplexId connectedEdgeId;
1163 do {
1164 // add a vertex
1165 const Cell vertex(0, currentId);
1166 vpath.push_back(vertex);
1167
1168 if(isCellCritical(vertex)) {
1169 break;
1170 }
1171
1172 connectedEdgeId = getPairedCell(vertex, triangulation);
1173 if(connectedEdgeId == -1) {
1174 break;
1175 }
1176
1177 // add an edge
1178 const Cell edge(1, connectedEdgeId);
1179 vpath.push_back(edge);
1180
1181 if(isCellCritical(edge)) {
1182 break;
1183 }
1184
1185 for(int i = 0; i < 2; ++i) {
1186 SimplexId vertexId;
1187 triangulation.getEdgeVertex(connectedEdgeId, i, vertexId);
1188
1189 if(vertexId != currentId) {
1190 currentId = vertexId;
1191 break;
1192 }
1193 }
1194
1195 } while(connectedEdgeId != -1);
1196 }
1197
1198 return 0;
1199}
1200
1201template <typename triangulationType>
1203 const Cell &saddle2,
1204 const Cell &saddle1,
1205 const std::vector<bool> &isVisited,
1206 std::vector<Cell> *const vpath,
1207 const triangulationType &triangulation,
1208 const bool stopIfMultiConnected,
1209 const bool enableCycleDetector) const {
1210
1211 // debug
1212 const SimplexId numberOfEdges = triangulation.getNumberOfEdges();
1213 std::vector<bool> isCycle;
1214 if(enableCycleDetector) {
1215 isCycle.resize(numberOfEdges, false);
1216 }
1217
1218 if(dimensionality_ == 3) {
1219 // add the 2-saddle to the path
1220 if(vpath != nullptr) {
1221 vpath->push_back(saddle2);
1222 }
1223
1224 SimplexId currentId = -1;
1225 {
1226 int nconnections = 0;
1227 for(int i = 0; i < 3; ++i) {
1228 SimplexId edgeId;
1229 triangulation.getTriangleEdge(saddle2.id_, i, edgeId);
1230 if(isVisited[edgeId]) {
1231 // saddle2 can be adjacent to saddle1 on the wall
1232 if(isSaddle1(Cell(1, edgeId))) {
1233 if(vpath != nullptr) {
1234 vpath->push_back(Cell(1, edgeId));
1235 }
1236 return false;
1237 }
1238
1239 currentId = edgeId;
1240 ++nconnections;
1241 }
1242 }
1243 if(stopIfMultiConnected && nconnections > 1) {
1244 return true;
1245 }
1246 }
1247
1248 int oldId;
1249 do {
1250
1251 // debug
1252 if(enableCycleDetector) {
1253 if(not isCycle[currentId]) {
1254 isCycle[currentId] = true;
1255 } else {
1256 this->printErr("Cycle detected on the wall of 1-saddle "
1257 + std::to_string(saddle1.id_));
1258 break;
1259 }
1260 }
1261
1262 oldId = currentId;
1263
1264 // add an edge
1265 const Cell edge(1, currentId);
1266 if(vpath != nullptr) {
1267 vpath->push_back(edge);
1268 }
1269
1270 if(isCellCritical(edge)) {
1271 break;
1272 }
1273
1274 const SimplexId connectedTriangleId = getPairedCell(edge, triangulation);
1275
1276 // add a triangle
1277 const Cell triangle(2, connectedTriangleId);
1278 if(vpath != nullptr) {
1279 vpath->push_back(triangle);
1280 }
1281
1282 if(isCellCritical(triangle)) {
1283 break;
1284 }
1285
1286 int nconnections = 0;
1287 for(int i = 0; i < 3; ++i) {
1288 SimplexId edgeId;
1289 triangulation.getTriangleEdge(connectedTriangleId, i, edgeId);
1290
1291 if(isVisited[edgeId] and edgeId != oldId) {
1292 currentId = edgeId;
1293 ++nconnections;
1294 }
1295 }
1296 if(stopIfMultiConnected && nconnections > 1) {
1297 return true;
1298 }
1299
1300 // stop at convergence caused by boundary effect
1301 } while(currentId != oldId);
1302 }
1303
1304 return false;
1305}
1306
1307template <typename triangulationType>
1309 std::vector<Cell> &vpath,
1310 const triangulationType &triangulation,
1311 const bool enableCycleDetector) const {
1312
1313 const SimplexId numberOfCells = triangulation.getNumberOfCells();
1314 std::vector<bool> isCycle;
1315 if(enableCycleDetector) {
1316 isCycle.resize(numberOfCells, false);
1317 }
1318
1319 if(dimensionality_ == 2) {
1320 if(cell.dim_ == 2) {
1321 // assume that cellId is a triangle
1322 SimplexId currentId = cell.id_;
1323 SimplexId oldId;
1324 do {
1325 oldId = currentId;
1326
1327 // add a triangle
1328 const Cell triangle(2, currentId);
1329 vpath.push_back(triangle);
1330
1331 if(isCellCritical(triangle)) {
1332 break;
1333 }
1334
1335 const SimplexId connectedEdgeId
1336 = getPairedCell(triangle, triangulation, true);
1337 if(connectedEdgeId == -1) {
1338 break;
1339 }
1340
1341 // add an edge
1342 const Cell edge(1, connectedEdgeId);
1343 vpath.push_back(edge);
1344
1345 if(isCellCritical(edge)) {
1346 break;
1347 }
1348
1349 const SimplexId starNumber
1350 = triangulation.getEdgeStarNumber(connectedEdgeId);
1351 for(SimplexId i = 0; i < starNumber; ++i) {
1352 SimplexId starId;
1353 triangulation.getEdgeStar(connectedEdgeId, i, starId);
1354
1355 if(starId != currentId) {
1356 currentId = starId;
1357 break;
1358 }
1359 }
1360
1361 // stop at convergence caused by boundary effect
1362 } while(currentId != oldId);
1363 }
1364 } else if(dimensionality_ == 3) {
1365 if(cell.dim_ == 3) {
1366 // assume that cellId is a tetra
1367 SimplexId currentId = cell.id_;
1368 SimplexId oldId;
1369 do {
1370
1371 // debug
1372 if(enableCycleDetector) {
1373 if(not isCycle[currentId]) {
1374 isCycle[currentId] = true;
1375 } else {
1376 this->printErr("cycle detected in the path from tetra "
1377 + std::to_string(cell.id_));
1378 break;
1379 }
1380 }
1381
1382 oldId = currentId;
1383
1384 // add a tetra
1385 const Cell tetra(3, currentId);
1386 vpath.push_back(tetra);
1387
1388 if(isCellCritical(tetra)) {
1389 break;
1390 }
1391
1392 const SimplexId connectedTriangleId
1393 = getPairedCell(tetra, triangulation, true);
1394 if(connectedTriangleId == -1) {
1395 break;
1396 }
1397
1398 // add a triangle
1399 const Cell triangle(2, connectedTriangleId);
1400 vpath.push_back(triangle);
1401
1402 if(isCellCritical(triangle)) {
1403 break;
1404 }
1405
1406 const SimplexId starNumber
1407 = triangulation.getTriangleStarNumber(connectedTriangleId);
1408 for(SimplexId i = 0; i < starNumber; ++i) {
1409 SimplexId starId;
1410 triangulation.getTriangleStar(connectedTriangleId, i, starId);
1411
1412 if(starId != currentId) {
1413 currentId = starId;
1414 break;
1415 }
1416 }
1417
1418 // stop at convergence caused by boundary effect
1419 } while(currentId != oldId);
1420 }
1421 }
1422
1423 return 0;
1424}
1425
1426template <typename triangulationType>
1428 const Cell &saddle1,
1429 const Cell &saddle2,
1430 const std::vector<bool> &isVisited,
1431 std::vector<Cell> *const vpath,
1432 const triangulationType &triangulation,
1433 const bool stopIfMultiConnected,
1434 const bool enableCycleDetector,
1435 bool *const cycleFound) const {
1436
1437 // debug
1438 const SimplexId numberOfTriangles = triangulation.getNumberOfTriangles();
1439 std::vector<bool> isCycle;
1440 if(enableCycleDetector) {
1441 isCycle.resize(numberOfTriangles, false);
1442 }
1443
1444 if(dimensionality_ == 3) {
1445 // add the 1-saddle to the path
1446 if(vpath != nullptr) {
1447 vpath->push_back(saddle1);
1448 }
1449
1450 SimplexId currentId = -1;
1451 {
1452 int nconnections = 0;
1453 const SimplexId triangleNumber
1454 = triangulation.getEdgeTriangleNumber(saddle1.id_);
1455 for(SimplexId i = 0; i < triangleNumber; ++i) {
1456 SimplexId triangleId;
1457 triangulation.getEdgeTriangle(saddle1.id_, i, triangleId);
1458 if(isVisited[triangleId]) {
1459 // saddle1 can be adjacent to saddle2 on the wall
1460 if(isSaddle2(Cell(2, triangleId))) {
1461 if(vpath != nullptr) {
1462 vpath->push_back(Cell(2, triangleId));
1463 }
1464 return false;
1465 }
1466
1467 currentId = triangleId;
1468 ++nconnections;
1469 }
1470 }
1471 if(stopIfMultiConnected && nconnections > 1) {
1472 return true;
1473 }
1474 }
1475
1476 if(currentId == -1) {
1477 return true;
1478 }
1479
1480 SimplexId oldId;
1481 do {
1482
1483 // debug
1484 if(enableCycleDetector) {
1485 if(not isCycle[currentId]) {
1486 isCycle[currentId] = true;
1487 } else {
1488 if(cycleFound)
1489 *cycleFound = true;
1490 else
1491 this->printErr("Cycle detected on the wall of 2-saddle "
1492 + std::to_string(saddle2.id_));
1493 break;
1494 }
1495 }
1496
1497 oldId = currentId;
1498
1499 // add a triangle
1500 const Cell triangle(2, currentId);
1501 if(vpath != nullptr) {
1502 vpath->push_back(triangle);
1503 }
1504
1505 if(isCellCritical(triangle)) {
1506 break;
1507 }
1508
1509 const SimplexId connectedEdgeId
1510 = getPairedCell(triangle, triangulation, true);
1511
1512 // add an edge
1513 const Cell edge(1, connectedEdgeId);
1514 if(vpath != nullptr) {
1515 vpath->push_back(edge);
1516 }
1517
1518 if(isCellCritical(edge)) {
1519 break;
1520 }
1521
1522 int nconnections = 0;
1523 const SimplexId triangleNumber
1524 = triangulation.getEdgeTriangleNumber(connectedEdgeId);
1525 for(SimplexId i = 0; i < triangleNumber; ++i) {
1526 SimplexId triangleId;
1527 triangulation.getEdgeTriangle(connectedEdgeId, i, triangleId);
1528
1529 if(isVisited[triangleId] and triangleId != oldId) {
1530 currentId = triangleId;
1531 ++nconnections;
1532 }
1533 }
1534 if(stopIfMultiConnected && nconnections > 1) {
1535 return true;
1536 }
1537
1538 // stop at convergence caused by boundary effect
1539 } while(currentId != oldId);
1540 }
1541
1542 return false;
1543}
1544
1545template <typename triangulationType>
1547 const Cell &cell, const triangulationType &triangulation) const {
1548 if(dimensionality_ == 3) {
1549 if(cell.dim_ == 1) {
1550 const SimplexId originId = getPairedCell(cell, triangulation);
1551 if(originId == -1)
1552 return false;
1553
1554 std::queue<SimplexId> bfs;
1555 bfs.push(originId);
1556 std::vector<bool> isVisited(triangulation.getNumberOfTriangles(), false);
1557
1558 // BFS traversal
1559 while(!bfs.empty()) {
1560 const SimplexId triangleId = bfs.front();
1561 bfs.pop();
1562
1563 isVisited[triangleId] = true;
1564
1565 for(int j = 0; j < 3; ++j) {
1566 SimplexId edgeId;
1567 triangulation.getTriangleEdge(triangleId, j, edgeId);
1568
1569 const SimplexId pairedCellId
1570 = getPairedCell(Cell(1, edgeId), triangulation);
1571
1572 if(triangleId == pairedCellId or pairedCellId == -1)
1573 continue;
1574
1575 if(isVisited[pairedCellId])
1576 return true;
1577 else
1578 bfs.push(pairedCellId);
1579 }
1580 }
1581 }
1582 }
1583 return false;
1584}
1585
1586template <typename triangulationType>
1588 const Cell &cell,
1589 VisitedMask &mask,
1590 const triangulationType &triangulation,
1591 std::vector<Cell> *const wall,
1592 std::vector<SimplexId> *const saddles) const {
1593
1594 if(saddles != nullptr) {
1595 saddles->clear();
1596 }
1597
1598 if(dimensionality_ == 3) {
1599 if(cell.dim_ == 2) {
1600 // assume that cellId is a triangle
1601 const SimplexId originId = cell.id_;
1602
1603 std::queue<SimplexId> bfs;
1604 bfs.push(originId);
1605
1606 // BFS traversal
1607 while(!bfs.empty()) {
1608 const SimplexId triangleId = bfs.front();
1609 bfs.pop();
1610
1611 if(!mask.isVisited_[triangleId]) {
1612 mask.isVisited_[triangleId] = true;
1613 mask.visitedIds_.emplace_back(triangleId);
1614
1615 // add the triangle
1616 if(wall != nullptr) {
1617 wall->push_back(Cell(2, triangleId));
1618 }
1619
1620 for(int j = 0; j < 3; ++j) {
1621 SimplexId edgeId;
1622 triangulation.getTriangleEdge(triangleId, j, edgeId);
1623
1624 if((saddles != nullptr) and isSaddle1(Cell(1, edgeId))) {
1625 saddles->emplace_back(edgeId);
1626 }
1627
1628 const SimplexId pairedCellId
1629 = getPairedCell(Cell(1, edgeId), triangulation);
1630
1631 if(pairedCellId != -1 and pairedCellId != triangleId) {
1632 bfs.push(pairedCellId);
1633 }
1634 }
1635 }
1636 }
1637
1638 if(saddles != nullptr && saddles->size() > 1) {
1639 std::sort(saddles->begin(), saddles->end());
1640 const auto last = std::unique(saddles->begin(), saddles->end());
1641 saddles->erase(last, saddles->end());
1642 }
1643 }
1644 }
1645
1646 return 0;
1647}
1648
1649template <typename triangulationType>
1651 const Cell &cell,
1652 VisitedMask &mask,
1653 const triangulationType &triangulation,
1654 std::vector<Cell> *const wall,
1655 std::vector<SimplexId> *const saddles) const {
1656
1657 if(saddles != nullptr) {
1658 saddles->clear();
1659 }
1660
1661 if(dimensionality_ == 3) {
1662 if(cell.dim_ == 1) {
1663 // assume that cellId is an edge
1664 const SimplexId originId = cell.id_;
1665
1666 std::queue<SimplexId> bfs;
1667 bfs.push(originId);
1668
1669 // BFS traversal
1670 while(!bfs.empty()) {
1671 const SimplexId edgeId = bfs.front();
1672 bfs.pop();
1673
1674 if(!mask.isVisited_[edgeId]) {
1675 mask.isVisited_[edgeId] = true;
1676 mask.visitedIds_.emplace_back(edgeId);
1677
1678 // add the edge
1679 if(wall != nullptr) {
1680 wall->push_back(Cell(1, edgeId));
1681 }
1682
1683 const SimplexId triangleNumber
1684 = triangulation.getEdgeTriangleNumber(edgeId);
1685 for(SimplexId j = 0; j < triangleNumber; ++j) {
1686 SimplexId triangleId;
1687 triangulation.getEdgeTriangle(edgeId, j, triangleId);
1688
1689 if((saddles != nullptr) and isSaddle2(Cell(2, triangleId))) {
1690 saddles->emplace_back(triangleId);
1691 }
1692
1693 const SimplexId pairedCellId
1694 = getPairedCell(Cell(2, triangleId), triangulation, true);
1695
1696 if(pairedCellId != -1 and pairedCellId != edgeId) {
1697 bfs.push(pairedCellId);
1698 }
1699 }
1700 }
1701 }
1702
1703 if(saddles != nullptr && saddles->size() > 1) {
1704 std::sort(saddles->begin(), saddles->end());
1705 const auto last = std::unique(saddles->begin(), saddles->end());
1706 saddles->erase(last, saddles->end());
1707 }
1708 }
1709 }
1710
1711 return 0;
1712}
1713
1714template <typename triangulationType>
1716 const std::vector<Cell> &vpath,
1717 const triangulationType &triangulation) const {
1718
1719 if(dimensionality_ == 2) {
1720 // assume that the first cell is an edge
1721 const SimplexId numberOfCellsInPath = vpath.size();
1722 for(SimplexId i = 0; i < numberOfCellsInPath; i += 2) {
1723 const SimplexId edgeId = vpath[i].id_;
1724 const SimplexId triangleId = vpath[i + 1].id_;
1725
1726#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1727 for(int k = 0; k < 3; ++k) {
1728 SimplexId tmp;
1729 triangulation.getCellEdge(triangleId, k, tmp);
1730 if(tmp == edgeId) {
1731 (*gradient_)[3][triangleId] = k;
1732 break;
1733 }
1734 }
1735 for(int k = 0; k < triangulation.getEdgeStarNumber(edgeId); ++k) {
1736 SimplexId tmp;
1737 triangulation.getEdgeStar(edgeId, k, tmp);
1738 if(tmp == triangleId) {
1739 (*gradient_)[2][edgeId] = k;
1740 break;
1741 }
1742 }
1743#else
1744 TTK_FORCE_USE(triangulation);
1745 (*gradient_)[3][triangleId] = edgeId;
1746 (*gradient_)[2][edgeId] = triangleId;
1747#endif
1748 }
1749 } else if(dimensionality_ == 3) {
1750 // assume that the first cell is a triangle
1751 const SimplexId numberOfCellsInPath = vpath.size();
1752 for(SimplexId i = 0; i < numberOfCellsInPath; i += 2) {
1753 const SimplexId triangleId = vpath[i].id_;
1754 const SimplexId tetraId = vpath[i + 1].id_;
1755
1756#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1757 for(int k = 0; k < 4; ++k) {
1758 SimplexId tmp;
1759 triangulation.getCellTriangle(tetraId, k, tmp);
1760 if(tmp == triangleId) {
1761 (*gradient_)[5][tetraId] = k;
1762 break;
1763 }
1764 }
1765 for(int k = 0; k < triangulation.getTriangleStarNumber(triangleId); ++k) {
1766 SimplexId tmp;
1767 triangulation.getTriangleStar(triangleId, k, tmp);
1768 if(tmp == tetraId) {
1769 (*gradient_)[4][triangleId] = k;
1770 break;
1771 }
1772 }
1773#else
1774 (*gradient_)[5][tetraId] = triangleId;
1775 (*gradient_)[4][triangleId] = tetraId;
1776#endif
1777 }
1778 }
1779
1780 return 0;
1781}
1782
1783template <typename triangulationType>
1785 const std::vector<Cell> &vpath,
1786 const triangulationType &triangulation) const {
1787
1788 // assume that the first cell is an edge
1789 for(size_t i = 0; i < vpath.size(); i += 2) {
1790 const SimplexId edgeId = vpath[i].id_;
1791 const SimplexId vertId = vpath[i + 1].id_;
1792
1793#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1794 const auto nneighs = triangulation.getVertexEdgeNumber();
1795 for(int k = 0; k < nneighs; ++k) {
1796 SimplexId tmp;
1797 triangulation.getVertexEdge(vertId, k, tmp);
1798 if(tmp == edgeId) {
1799 (*gradient_)[0][vertId] = k;
1800 break;
1801 }
1802 }
1803 const auto nverts = triangulation.getEdgeStarNumber(edgeId);
1804 for(int k = 0; k < nverts; ++k) {
1805 SimplexId tmp;
1806 triangulation.getEdgeVertex(edgeId, k, tmp);
1807 if(tmp == vertId) {
1808 (*gradient_)[1][edgeId] = k;
1809 break;
1810 }
1811 }
1812#else
1813 TTK_FORCE_USE(triangulation);
1814 (*gradient_)[0][vertId] = edgeId;
1815 (*gradient_)[1][edgeId] = vertId;
1816#endif
1817 }
1818
1819 return 0;
1820}
1821
1822template <typename triangulationType>
1824 const std::vector<Cell> &vpath,
1825 const triangulationType &triangulation,
1826 bool cancelReversal) const {
1827
1828 if(dimensionality_ == 3) {
1829 // assume that the first cell is an edge
1830 if(cancelReversal) {
1831 if(vpath.empty())
1832 return 0;
1833 (*gradient_)[2][vpath[0].id_] = NULL_GRADIENT;
1834 // assume that the last cell is a triangle
1835 if(vpath.size() <= 1)
1836 return 0;
1837 (*gradient_)[3][vpath[vpath.size() - 1].id_] = NULL_GRADIENT;
1838 }
1839 const SimplexId numberOfCellsInPath = vpath.size();
1840 const SimplexId startIndex = (cancelReversal ? 2 : 0);
1841 for(SimplexId i = startIndex; i < numberOfCellsInPath; i += 2) {
1842 const SimplexId vpathEdgeIndex = i;
1843 const SimplexId vpathTriangleIndex = (cancelReversal ? i - 1 : i + 1);
1844 const SimplexId edgeId = vpath[vpathEdgeIndex].id_;
1845 const SimplexId triangleId = vpath[vpathTriangleIndex].id_;
1846
1847#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1848 for(int k = 0; k < 3; ++k) {
1849 SimplexId tmp;
1850 triangulation.getTriangleEdge(triangleId, k, tmp);
1851 if(tmp == edgeId) {
1852 (*gradient_)[3][triangleId] = k;
1853 break;
1854 }
1855 }
1856 for(int k = 0; k < triangulation.getEdgeTriangleNumber(edgeId); ++k) {
1857 SimplexId tmp;
1858 triangulation.getEdgeTriangle(edgeId, k, tmp);
1859 if(tmp == triangleId) {
1860 (*gradient_)[2][edgeId] = k;
1861 break;
1862 }
1863 }
1864#else
1865 TTK_FORCE_USE(triangulation);
1866 (*gradient_)[3][triangleId] = edgeId;
1867 (*gradient_)[2][edgeId] = triangleId;
1868#endif
1869 }
1870 }
1871
1872 return 0;
1873}
1874
1875template <typename triangulationType>
1877 const std::vector<Cell> &vpath,
1878 const triangulationType &triangulation) const {
1879
1880 if(dimensionality_ == 3) {
1881 // assume that the first cell is a triangle
1882 const SimplexId numberOfCellsInPath = vpath.size();
1883 for(SimplexId i = 0; i < numberOfCellsInPath; i += 2) {
1884 const SimplexId triangleId = vpath[i].id_;
1885 const SimplexId edgeId = vpath[i + 1].id_;
1886
1887#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
1888 for(int k = 0; k < 3; ++k) {
1889 SimplexId tmp;
1890 triangulation.getTriangleEdge(triangleId, k, tmp);
1891 if(tmp == edgeId) {
1892 (*gradient_)[3][triangleId] = k;
1893 break;
1894 }
1895 }
1896 for(int k = 0; k < triangulation.getEdgeTriangleNumber(edgeId); ++k) {
1897 SimplexId tmp;
1898 triangulation.getEdgeTriangle(edgeId, k, tmp);
1899 if(tmp == triangleId) {
1900 (*gradient_)[2][edgeId] = k;
1901 break;
1902 }
1903 }
1904#else
1905 TTK_FORCE_USE(triangulation);
1906 (*gradient_)[2][edgeId] = triangleId;
1907 (*gradient_)[3][triangleId] = edgeId;
1908#endif
1909 }
1910 }
1911
1912 return 0;
1913}
1914
1915template <typename triangulationType>
1917 const Cell c, const triangulationType &triangulation) const {
1918
1919 const auto offsets = this->inputOffsets_;
1920
1921 auto cellDim = c.dim_;
1922 auto cellId = c.id_;
1923
1924 SimplexId vertexId = -1;
1925 if(cellDim == 0) {
1926 vertexId = cellId;
1927 }
1928
1929 else if(cellDim == 1) {
1930 SimplexId v0;
1931 SimplexId v1;
1932 triangulation.getEdgeVertex(cellId, 0, v0);
1933 triangulation.getEdgeVertex(cellId, 1, v1);
1934
1935 if(offsets[v0] > offsets[v1]) {
1936 vertexId = v0;
1937 } else {
1938 vertexId = v1;
1939 }
1940 }
1941
1942 else if(cellDim == 2) {
1943 SimplexId v0{}, v1{}, v2{};
1944 triangulation.getTriangleVertex(cellId, 0, v0);
1945 triangulation.getTriangleVertex(cellId, 1, v1);
1946 triangulation.getTriangleVertex(cellId, 2, v2);
1947 if(offsets[v0] > offsets[v1] && offsets[v0] > offsets[v2]) {
1948 vertexId = v0;
1949 } else if(offsets[v1] > offsets[v0] && offsets[v1] > offsets[v2]) {
1950 vertexId = v1;
1951 } else {
1952 vertexId = v2;
1953 }
1954 }
1955
1956 else if(cellDim == 3) {
1957 SimplexId v0{}, v1{}, v2{}, v3{};
1958 triangulation.getCellVertex(cellId, 0, v0);
1959 triangulation.getCellVertex(cellId, 1, v1);
1960 triangulation.getCellVertex(cellId, 2, v2);
1961 triangulation.getCellVertex(cellId, 3, v3);
1962 if(offsets[v0] > offsets[v1] && offsets[v0] > offsets[v2]
1963 && offsets[v0] > offsets[v3]) {
1964 vertexId = v0;
1965 } else if(offsets[v1] > offsets[v0] && offsets[v1] > offsets[v2]
1966 && offsets[v1] > offsets[v3]) {
1967 vertexId = v1;
1968 } else if(offsets[v2] > offsets[v0] && offsets[v2] > offsets[v1]
1969 && offsets[v2] > offsets[v3]) {
1970 vertexId = v2;
1971 } else {
1972 vertexId = v3;
1973 }
1974 }
1975 return vertexId;
1976}
1977
1978template <typename triangulationType>
1980 const Cell c, const triangulationType &triangulation) const {
1981
1982 const auto offsets = this->inputOffsets_;
1983
1984 auto cellDim = c.dim_;
1985 auto cellId = c.id_;
1986
1987 SimplexId vertexId = -1;
1988 if(cellDim == 0) {
1989 vertexId = cellId;
1990 }
1991
1992 else if(cellDim == 1) {
1993 SimplexId v0;
1994 SimplexId v1;
1995 triangulation.getEdgeVertex(cellId, 0, v0);
1996 triangulation.getEdgeVertex(cellId, 1, v1);
1997
1998 if(offsets[v0] < offsets[v1]) {
1999 vertexId = v0;
2000 } else {
2001 vertexId = v1;
2002 }
2003 }
2004
2005 else if(cellDim == 2) {
2006 SimplexId v0{}, v1{}, v2{};
2007 triangulation.getTriangleVertex(cellId, 0, v0);
2008 triangulation.getTriangleVertex(cellId, 1, v1);
2009 triangulation.getTriangleVertex(cellId, 2, v2);
2010 if(offsets[v0] < offsets[v1] && offsets[v0] < offsets[v2]) {
2011 vertexId = v0;
2012 } else if(offsets[v1] < offsets[v0] && offsets[v1] < offsets[v2]) {
2013 vertexId = v1;
2014 } else {
2015 vertexId = v2;
2016 }
2017 }
2018
2019 else if(cellDim == 3) {
2020 SimplexId v0{}, v1{}, v2{}, v3{};
2021 triangulation.getCellVertex(cellId, 0, v0);
2022 triangulation.getCellVertex(cellId, 1, v1);
2023 triangulation.getCellVertex(cellId, 2, v2);
2024 triangulation.getCellVertex(cellId, 3, v3);
2025 if(offsets[v0] < offsets[v1] && offsets[v0] < offsets[v2]
2026 && offsets[v0] < offsets[v3]) {
2027 vertexId = v0;
2028 } else if(offsets[v1] < offsets[v0] && offsets[v1] < offsets[v2]
2029 && offsets[v1] < offsets[v3]) {
2030 vertexId = v1;
2031 } else if(offsets[v2] < offsets[v0] && offsets[v2] < offsets[v1]
2032 && offsets[v2] < offsets[v3]) {
2033 vertexId = v2;
2034 } else {
2035 vertexId = v3;
2036 }
2037 }
2038 return vertexId;
2039}
2040
2041template <typename triangulationType>
2043 std::vector<std::array<float, 3>> &points,
2044 std::vector<char> &points_pairOrigins,
2045 std::vector<char> &cells_pairTypes,
2046 std::vector<SimplexId> &cellIds,
2047 std::vector<char> &cellDimensions,
2048 const triangulationType &triangulation) const {
2049
2050 const auto nDims = this->getNumberOfDimensions();
2051
2052 // number of glyphs per dimension
2053 std::vector<size_t> nGlyphsPerDim(nDims);
2054
2055#ifdef TTK_ENABLE_OPENMP
2056#pragma omp parallel for num_threads(threadNumber_)
2057#endif // TTK_ENABLE_OPENMP
2058 for(int i = 0; i < nDims - 1; ++i) {
2059 const auto nCells = this->getNumberOfCells(i, triangulation);
2060 for(SimplexId j = 0; j < nCells; ++j) {
2061 if(this->getPairedCell(Cell{i, j}, triangulation) > -1) {
2062 nGlyphsPerDim[i]++;
2063 }
2064 }
2065 }
2066
2067 // partial sum of number of gradient glyphs
2068 std::vector<size_t> offsets(nDims + 1);
2069 for(SimplexId i = 0; i < nDims; ++i) {
2070 offsets[i + 1] = offsets[i] + nGlyphsPerDim[i];
2071 }
2072
2073 // total number of glyphs
2074 const auto nGlyphs = offsets.back();
2075
2076 // resize arrays accordingly
2077 points.resize(2 * nGlyphs);
2078 points_pairOrigins.resize(2 * nGlyphs);
2079 cells_pairTypes.resize(nGlyphs);
2080 cellIds.resize(2 * nGlyphs);
2081 cellDimensions.resize(2 * nGlyphs);
2082
2083#ifdef TTK_ENABLE_OPENMP
2084#pragma omp parallel for num_threads(threadNumber_)
2085#endif // TTK_ENABLE_OPENMP
2086 for(int i = 0; i < nDims - 1; ++i) {
2087 const SimplexId nCells = getNumberOfCells(i, triangulation);
2088 size_t nProcessedGlyphs{offsets[i]};
2089 for(SimplexId j = 0; j < nCells; ++j) {
2090 const Cell c{i, j};
2091 const auto pcid = this->getPairedCell(c, triangulation);
2092 if(pcid > -1) {
2093 const Cell pc{i + 1, pcid};
2094 triangulation.getCellIncenter(
2095 c.id_, c.dim_, points[2 * nProcessedGlyphs].data());
2096 triangulation.getCellIncenter(
2097 pc.id_, pc.dim_, points[2 * nProcessedGlyphs + 1].data());
2098 points_pairOrigins[2 * nProcessedGlyphs] = 0;
2099 points_pairOrigins[2 * nProcessedGlyphs + 1] = 1;
2100 cells_pairTypes[nProcessedGlyphs] = i;
2101#ifdef TTK_ENABLE_MPI
2102 ttk::SimplexId globalId{-1};
2103 triangulation.getDistributedGlobalCellId(j, i, globalId);
2104 cellIds[2 * nProcessedGlyphs + 0] = globalId;
2105 triangulation.getDistributedGlobalCellId(pcid, i + 1, globalId);
2106 cellIds[2 * nProcessedGlyphs + 1] = globalId;
2107#else
2108 cellIds[2 * nProcessedGlyphs + 0] = j;
2109 cellIds[2 * nProcessedGlyphs + 1] = pcid;
2110#endif // TTK_ENABLE_MPI
2111 cellDimensions[2 * nProcessedGlyphs + 0] = i;
2112 cellDimensions[2 * nProcessedGlyphs + 1] = i + 1;
2113 nProcessedGlyphs++;
2114 }
2115 }
2116 }
2117
2118 return 0;
2119}
#define TTK_FORCE_USE(x)
Force the compiler to use the function/method parameter.
Definition BaseClass.h:57
std::array< std::vector< gradIdType >, 6 > gradientType
Discrete gradient internal struct.
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
TTK discreteGradient processing package.
SimplexId getCellLowerVertex(const Cell c, const triangulationType &triangulation) const
int reverseAscendingPath(const std::vector< Cell > &vpath, const triangulationType &triangulation) const
int reverseDescendingPath(const std::vector< Cell > &vpath, const triangulationType &triangulation) const
int setGradientGlyphs(std::vector< std::array< float, 3 > > &points, std::vector< char > &points_pairOrigins, std::vector< char > &cells_pairTypes, std::vector< SimplexId > &cellsIds, std::vector< char > &cellsDimensions, const triangulationType &triangulation) const
int reverseAscendingPathOnWall(const std::vector< Cell > &vpath, const triangulationType &triangulation, bool cancelReversal=false) const
bool getDescendingPathThroughWall(const Cell &saddle2, const Cell &saddle1, const std::vector< bool > &isVisited, std::vector< Cell > *const vpath, const triangulationType &triangulation, const bool stopIfMultiConnected=false, const bool enableCycleDetector=false) const
bool isSaddle1(const Cell &cell) const
bool isSaddle2(const Cell &cell) const
AbstractTriangulation::gradientKeyType inputScalarField_
SimplexId getNumberOfCells(const int dimension, const triangulationType &triangulation) const
bool isBoundary(const Cell &cell, const triangulationType &triangulation) const
bool isCellCritical(const int cellDim, const SimplexId cellId) const
SimplexId getPairedCell(const Cell &cell, const triangulationType &triangulation, bool isReverse=false) const
int getAscendingWall(const Cell &cell, VisitedMask &mask, const triangulationType &triangulation, std::vector< Cell > *const wall=nullptr, std::vector< SimplexId > *const saddles=nullptr) const
int reverseDescendingPathOnWall(const std::vector< Cell > &vpath, const triangulationType &triangulation) const
SimplexId getCellGreaterVertex(const Cell c, const triangulationType &triangulation) const
int getCriticalPoints(std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, const triangulationType &triangulation) const
int buildGradient(const triangulationType &triangulation, bool bypassCache=false, const std::vector< bool > *updateMask=nullptr)
bool getAscendingPathThroughWall(const Cell &saddle1, const Cell &saddle2, const std::vector< bool > &isVisited, std::vector< Cell > *const vpath, const triangulationType &triangulation, const bool stopIfMultiConnected=false, const bool enableCycleDetector=false, bool *const cycleFound=nullptr) const
dataType getPersistence(const Cell &up, const Cell &down, const dataType *const scalars, const triangulationType &triangulation) const
int setCriticalPoints(const std::array< std::vector< SimplexId >, 4 > &criticalCellsByDim, std::vector< std::array< float, 3 > > &points, std::vector< char > &cellDimensions, std::vector< SimplexId > &cellIds, std::vector< char > &isOnBoundary, std::vector< SimplexId > &PLVertexIdentifiers, const triangulationType &triangulation) const
int getDescendingPath(const Cell &cell, std::vector< Cell > &vpath, const triangulationType &triangulation) const
bool detectGradientCycle(const Cell &cell, const triangulationType &triangulation) const
int getAscendingPath(const Cell &cell, std::vector< Cell > &vpath, const triangulationType &triangulation, const bool enableCycleDetector=false) const
AbstractTriangulation::gradientType localGradient_
int getDescendingWall(const Cell &cell, VisitedMask &mask, const triangulationType &triangulation, std::vector< Cell > *const wall=nullptr, std::vector< SimplexId > *const saddles=nullptr) const
AbstractTriangulation::gradientType * gradient_
const SimplexId * inputOffsets_
COMMON_EXPORTS int MPIrank_
Definition BaseClass.cpp:9
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
T end(std::pair< T, T > &p)
Definition ripserpy.cpp:483
Auto-cleaning re-usable graph propagations data structure.
Definition VisitedMask.h:27
std::vector< bool > & isVisited_
Definition VisitedMask.h:28
std::vector< SimplexId > & visitedIds_
Definition VisitedMask.h:29
Extended Cell structure for processLowerStars.
const std::array< SimplexId, 3 > lowVerts_
const std::array< uint8_t, 3 > faces_
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)