172 template <
typename dataType,
typename triangulationType>
177 const dataType *
const scalars,
178 const size_t scalarsMTime,
180 const triangulationType &triangulation,
181 const unsigned int &seed = 0);
195 const bool doDescending,
196 const bool doSaddleConnectors) {
206 const bool doDescending) {
215 const bool doDescending,
216 const bool doMorseSmale) {
283 template <
typename triangulationType>
286 std::vector<Separatrix> &separatrices,
287 const triangulationType &triangulation)
const;
293 template <
typename triangulationType>
295 std::vector<Separatrix> &separatrices,
296 const triangulationType &triangulation)
const;
302 template <
typename triangulationType>
304 std::vector<Separatrix> &separatrices,
305 const triangulationType &triangulation)
const;
310 template <
typename triangulationType>
312 const std::vector<Separatrix> &separatrices,
314 const triangulationType &triangulation)
const;
320 template <
typename triangulationType>
322 const std::vector<SimplexId> &saddles2,
323 std::vector<Separatrix> &separatrices,
324 std::vector<std::vector<SimplexId>> &separatricesSaddles,
325 const triangulationType &triangulation)
const;
331 template <
typename triangulationType>
334 const std::vector<Separatrix> &separatrices,
335 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
337 const triangulationType &triangulation)
const;
344 template <
typename triangulationType>
347 const size_t polSize,
348 const triangulationType &triangulation)
const;
353 template <
typename triangulationType>
355 const size_t polSize,
356 const triangulationType &triangulation)
const;
362 template <
typename triangulationType>
364 const std::vector<SimplexId> &saddles1,
365 std::vector<Separatrix> &separatrices,
366 std::vector<std::vector<SimplexId>> &separatricesSaddles,
367 const triangulationType &triangulation)
const;
373 template <
typename triangulationType>
376 const std::vector<Separatrix> &separatrices,
377 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
379 const triangulationType &triangulation)
const;
385 std::vector<std::vector<Separatrix>> &separatrices)
const;
390 template <
typename triangulationType>
393 const triangulationType &triangulation)
const;
398 template <
typename triangulationType>
401 const triangulationType &triangulation)
const;
407 template <
typename triangulationType>
409 const SimplexId *
const ascendingManifold,
410 const SimplexId *
const descendingManifold,
412 const triangulationType &triangulation)
const;
414 template <
typename dataType,
typename triangulationType>
416 const dataType *
const scalars,
418 const triangulationType &triangulation);
420 template <
typename dataType,
typename triangulationType>
422 const double &persistenceThreshold,
423 const bool &detectCriticalPoints
455 const dataType *
const scalars,
456 const size_t scalarsMTime,
458 const triangulationType &triangulation,
459 const unsigned int &seed) {
460#ifndef TTK_ENABLE_KAMIKAZE
461 if(scalars ==
nullptr) {
462 this->
printErr(
"Input scalar field pointer is null.");
466 if(offsets ==
nullptr) {
467 this->
printErr(
"Input offset field pointer is null.");
476 const auto dim = triangulation.getDimensionality();
484 == DiscreteGradient::BACKEND::STOCHASTIC_BACKEND) {
492 const auto nVerts{triangulation.getNumberOfVertices()};
495 const auto pair{std::minmax_element(offsets, offsets + nVerts)};
497 const auto globmin = std::distance(offsets, pair.first);
498 const auto globmax = std::distance(offsets, pair.second);
499 persistenceThreshold *= (scalars[globmax] - scalars[globmin]);
500 this->
printMsg(
"Absolute saddle connectors persistence threshold is "
501 + std::to_string(persistenceThreshold),
506 == DiscreteGradient::BACKEND::CLASSIC_BACKEND)
508 persistenceThreshold, scalars, offsets, triangulation);
510 == DiscreteGradient::BACKEND::STOCHASTIC_BACKEND) {
512 triangulation, persistenceThreshold,
true);
516 std::array<std::vector<SimplexId>, 4> criticalPoints{};
520 this->
printMsg(
" Critical points extracted", 1.0, tm.getElapsedTime(),
525 std::vector<std::vector<Separatrix>> separatrices1{};
532 separatrices1.emplace_back();
535 criticalPoints[1], separatrices1.back(), triangulation);
537 this->
printMsg(
" Descending 1-separatrices computed", 1.0,
544 separatrices1.emplace_back();
547 criticalPoints[dim - 1], separatrices1.back(), triangulation);
549 this->
printMsg(
" Ascending 1-separatrices computed", 1.0,
557 separatrices1.emplace_back();
574 this->
printMsg(
" 1-separatrices set", 1.0, tmp.getElapsedTime(),
578 this->
printMsg(
"1-separatrices computed", 1.0, tm1sep.getElapsedTime(),
579 this->threadNumber_);
587 std::vector<Separatrix> separatrices;
588 std::vector<std::vector<SimplexId>> separatricesSaddles;
590 criticalPoints[2], separatrices, separatricesSaddles, triangulation);
592 outSeps2, separatrices, separatricesSaddles, offsets, triangulation);
594 this->
printMsg(
" Descending 2-separatrices computed", 1.0,
601 std::vector<Separatrix> separatrices;
602 std::vector<std::vector<SimplexId>> separatricesSaddles;
604 criticalPoints[1], separatrices, separatricesSaddles, triangulation);
606 outSeps2, separatrices, separatricesSaddles, offsets, triangulation);
608 this->
printMsg(
" Ascending 2-separatrices computed", 1.0,
615 this->
printMsg(
"2-separatrices computed", 1.0, tm2sep.getElapsedTime(),
616 this->threadNumber_);
624 criticalPoints[dim], outManifold.
ascending_, triangulation);
628 criticalPoints[0], outManifold.
descending_, triangulation);
638 "Segmentation computed", 1.0, tmp.
getElapsedTime(), this->threadNumber_);
654 + std::to_string(triangulation.getNumberOfVertices())
655 +
" points) processed",
712 const std::vector<SimplexId> &saddles,
713 std::vector<Separatrix> &separatrices,
714 const triangulationType &triangulation)
const {
719 const auto dim{triangulation.getDimensionality()};
722 auto getFaceStarNumber = &triangulationType::getTriangleStarNumber;
723 auto getFaceStar = &triangulationType::getTriangleStar;
726 getFaceStarNumber = &triangulationType::getEdgeStarNumber;
727 getFaceStar = &triangulationType::getEdgeStar;
730 const SimplexId numberOfSaddles = saddles.size();
732 std::vector<std::vector<Separatrix>> sepsPerSaddle(numberOfSaddles);
735#ifdef TTK_ENABLE_OPENMP
736#pragma omp parallel for num_threads(threadNumber_)
738 for(
SimplexId i = 0; i < numberOfSaddles; ++i) {
739 const Cell saddle{dim - 1, saddles[i]};
742 const auto starNumber{(triangulation.*getFaceStarNumber)(saddle.id_)};
743 for(
SimplexId j = 0; j < starNumber; ++j) {
746 (triangulation.*getFaceStar)(saddle.id_, j, sId);
748 std::vector<Cell> vpath{saddle};
751 const Cell &lastCell = vpath.back();
753 sepsPerSaddle[i].emplace_back();
754 sepsPerSaddle[i].back().source_ = saddle;
755 sepsPerSaddle[i].back().destination_ = lastCell;
756 sepsPerSaddle[i].back().geometry_ = std::move(vpath);
763 separatrices = std::move(sepsPerSaddle[0]);
770 const std::vector<SimplexId> &saddles2,
771 std::vector<Separatrix> &separatrices,
772 const triangulationType &triangulation)
const {
777 const auto nTriangles = triangulation.getNumberOfTriangles();
779 std::vector<bool> isVisited(nTriangles,
false);
780 std::vector<SimplexId> visitedTriangles{};
782 using Vpath = std::vector<Cell>;
784 const auto dim{triangulation.getDimensionality()};
786 std::vector<std::vector<Separatrix>> sepsByThread(saddles2.size());
787 std::vector<SimplexId> saddles1{};
789#ifdef TTK_ENABLE_OPENMP
790#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
791 firstprivate(isVisited, visitedTriangles, saddles1)
793 for(
size_t i = 0; i < saddles2.size(); ++i) {
794 const Cell s2{dim - 1, saddles2[i]};
798 s2, mask, triangulation,
nullptr, &saddles1);
800 for(
const auto saddle1Id : saddles1) {
801 const Cell s1{1, saddle1Id};
804 const bool isMultiConnected
806 s1, s2, isVisited, &vpath, triangulation);
812 const auto &last = vpath.back();
814 if(!isMultiConnected && last.dim_ == s2.dim_ && last.id_ == s2.id_) {
815 sepsByThread[i].emplace_back();
816 sepsByThread[i].back().source_ = s1;
817 sepsByThread[i].back().destination_ = s2;
818 sepsByThread[i].back().geometry_ = std::move(vpath);
825 separatrices = std::move(sepsByThread[0]);
833 const std::vector<Separatrix> &separatrices,
835 const triangulationType &triangulation)
const {
853 std::vector<size_t> geomPointsBegId{npoints};
855 std::vector<size_t> geomCellsBegId{ncells};
858 for(
size_t i = 0; i < separatrices.size(); ++i) {
859 const auto &sep = separatrices[i];
860 const auto sepSize = sep.geometry_.size();
862 ncells += sepSize - 1;
863 geomPointsBegId.emplace_back(npoints);
864 geomCellsBegId.emplace_back(ncells);
867 const int dimensionality = triangulation.getCellVertexNumber(0) - 1;
881 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
882 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
885#ifdef TTK_ENABLE_OPENMP
886#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic)
888 for(
size_t i = 0; i < separatrices.size(); ++i) {
889 const auto &sep = separatrices[i];
890 const auto &sepGeom = sep.geometry_;
891 const auto sepId = separatrixId + i;
898 const auto saddleConnector
899 = dimensionality == 3 && src.
dim_ == 1 && dst.
dim_ == 2;
901 = saddleConnector ? 1 : std::min(dst.
dim_, dimensionality - 1);
905 return offsets[a] < offsets[b];
907 const std::array<SimplexId, 2> gVerts{
910 const auto sepFuncMax
911 = *std::max_element(gVerts.begin(), gVerts.end(), vertsOrder);
912 const std::array<SimplexId, 2> lVerts{
915 const auto sepFuncMin
916 = *std::min_element(lVerts.begin(), lVerts.end(), vertsOrder);
917 separatrixFunctionMaxima[sepId] = sepFuncMax;
918 separatrixFunctionMinima[sepId] = sepFuncMin;
921 const auto onBoundary
925 for(
size_t j = 0; j < sepGeom.size(); ++j) {
926 const auto &cell = sepGeom[j];
927 std::array<float, 3> pt{};
928 triangulation.getCellIncenter(cell.id_, cell.dim_, pt.data());
931 const auto k = geomPointsBegId[i] + j;
933 points[3 * k + 0] = pt[0];
934 points[3 * k + 1] = pt[1];
935 points[3 * k + 2] = pt[2];
938 = (j == 0 || j == sepGeom.size() - 1) ? 0 : 1;
947 const auto l = geomCellsBegId[i] + j - 1;
949 cellsConn[2 * l + 0] = k - 1;
950 cellsConn[2 * l + 1] = k;
973 const std::vector<SimplexId> &saddles1,
974 std::vector<Separatrix> &separatrices,
975 std::vector<std::vector<SimplexId>> &separatricesSaddles,
976 const triangulationType &triangulation)
const {
977 const Cell emptyCell;
979 const SimplexId numberOfSaddles = saddles1.size();
983 const SimplexId numberOfSeparatrices = numberOfSaddles;
984 separatrices.resize(numberOfSeparatrices);
985 separatricesSaddles.resize(numberOfSeparatrices);
987 const auto nEdges = triangulation.getNumberOfEdges();
988 std::vector<bool> isVisited(nEdges,
false);
989 std::vector<SimplexId> visitedEdges{};
992#ifdef TTK_ENABLE_OPENMP
993#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
994 firstprivate(isVisited, visitedEdges)
996 for(
SimplexId i = 0; i < numberOfSaddles; ++i) {
997 const Cell saddle1{1, saddles1[i]};
999 std::vector<Cell> wall;
1002 saddle1, mask, triangulation, &wall, &separatricesSaddles[i]);
1004 separatrices[i].source_ = saddle1;
1005 separatrices[i].destination_ = emptyCell;
1006 separatrices[i].geometry_ = std::move(wall);
1014 const std::vector<SimplexId> &saddles2,
1015 std::vector<Separatrix> &separatrices,
1016 std::vector<std::vector<SimplexId>> &separatricesSaddles,
1017 const triangulationType &triangulation)
const {
1018 const Cell emptyCell;
1020 const SimplexId numberOfSaddles = saddles2.size();
1024 const SimplexId numberOfSeparatrices = numberOfSaddles;
1025 separatrices.resize(numberOfSeparatrices);
1026 separatricesSaddles.resize(numberOfSeparatrices);
1028 const auto nTriangles = triangulation.getNumberOfTriangles();
1029 std::vector<bool> isVisited(nTriangles,
false);
1030 std::vector<SimplexId> visitedTriangles{};
1032 const auto dim{triangulation.getDimensionality()};
1035#ifdef TTK_ENABLE_OPENMP
1036#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic) \
1037 firstprivate(isVisited, visitedTriangles)
1039 for(
SimplexId i = 0; i < numberOfSaddles; ++i) {
1040 const Cell saddle2{dim - 1, saddles2[i]};
1042 std::vector<Cell> wall;
1045 saddle2, mask, triangulation, &wall, &separatricesSaddles[i]);
1047 separatrices[i].source_ = saddle2;
1048 separatrices[i].destination_ = emptyCell;
1049 separatrices[i].geometry_ = std::move(wall);
1109 const std::vector<Separatrix> &separatrices,
1110 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
1112 const triangulationType &triangulation)
const {
1130 const auto noldcells{ncells};
1134 std::vector<size_t> geomCellsBegId{ncells};
1137 for(
size_t i = 0; i < separatrices.size(); ++i) {
1138 const auto &sep = separatrices[i];
1139 ncells += sep.geometry_.size();
1140 geomCellsBegId.emplace_back(ncells);
1144 std::vector<SimplexId> sepSourceIds(separatrices.size());
1145 std::vector<SimplexId> sepIds(separatrices.size());
1146 std::vector<char> sepOnBoundary(separatrices.size());
1147 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
1148 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
1150 std::vector<SimplexId> polygonNTetras(ncells - noldcells);
1151 std::vector<SimplexId> polygonEdgeIds(ncells - noldcells);
1152 std::vector<SimplexId> polygonSepInfosIds(ncells - noldcells);
1154#ifdef TTK_ENABLE_OPENMP
1155#pragma omp parallel for num_threads(threadNumber_) schedule(dynamic)
1157 for(
size_t i = 0; i < separatrices.size(); ++i) {
1158 const auto &sep = separatrices[i];
1159 const auto &sepGeom = sep.geometry_;
1160 const auto &sepSaddles = separatricesSaddles[i];
1161 const auto sepId = separatrixId + i;
1165 const auto sepFuncMin
1168 if(!sepSaddles.empty()) {
1170 const auto maxId = *std::max_element(
1171 sepSaddles.begin(), sepSaddles.end(),
1173 return offsets[discreteGradient_.getCellGreaterVertex(
1174 Cell{2, a}, triangulation)]
1176 Cell{2, b}, triangulation)];
1183 const auto maxId = *std::max_element(
1184 sepGeom.begin(), sepGeom.end(),
1185 [&triangulation, offsets,
this](
const Cell &a,
const Cell &b) {
1186 return offsets[discreteGradient_.getCellGreaterVertex(
1188 < offsets[discreteGradient_.getCellGreaterVertex(
1195 const char onBoundary
1196 = (sepSaddles.empty()
1198 : std::count_if(sepSaddles.begin(), sepSaddles.end(),
1200 return triangulation.isTriangleOnBoundary(a);
1202 + triangulation.isEdgeOnBoundary(src.id_);
1205 sepSourceIds[i] = src.id_;
1206 separatrixFunctionMaxima[sepId] = sepFuncMax;
1207 separatrixFunctionMinima[sepId] = sepFuncMin;
1208 sepOnBoundary[i] = onBoundary;
1210 for(
size_t j = 0; j < sepGeom.size(); ++j) {
1211 const auto &cell = sepGeom[j];
1213 const auto k = geomCellsBegId[i] + j - noldcells;
1215 polygonNTetras[k] = triangulation.getEdgeStarNumber(cell.id_);
1217 if(polygonNTetras[k] > 2) {
1218 polygonEdgeIds[k] = cell.id_;
1219 polygonSepInfosIds[k] = i;
1225 std::vector<SimplexId> validTetraIds{};
1226 validTetraIds.reserve(polygonNTetras.size());
1228 for(
size_t i = 0; i < polygonNTetras.size(); ++i) {
1229 if(polygonNTetras[i] > 2) {
1230 validTetraIds.emplace_back(i);
1235 size_t nnewpoints{};
1236 std::vector<SimplexId> pointsPerCell(validTetraIds.size() + 1);
1237 for(
size_t i = 0; i < validTetraIds.size(); ++i) {
1238 nnewpoints += polygonNTetras[validTetraIds[i]];
1239 pointsPerCell[i + 1] = nnewpoints;
1243 outSeps2.cl.connectivity_.resize(firstCellId + nnewpoints);
1244 auto cellsConn = &outSeps2.cl.connectivity_[firstCellId];
1246 std::vector<SimplexId> cellVertsIds(nnewpoints);
1248#ifdef TTK_ENABLE_OPENMP
1249#pragma omp parallel for num_threads(threadNumber_)
1251 for(
size_t i = 0; i < validTetraIds.size(); ++i) {
1252 const auto k = validTetraIds[i];
1255 getDualPolygon(polygonEdgeIds[k], &cellVertsIds[pointsPerCell[i]],
1256 polygonNTetras[k], triangulation);
1258 sortDualPolygonVertices(
1259 &cellVertsIds[pointsPerCell[i]], polygonNTetras[k], triangulation);
1261 for(
SimplexId j = 0; j < polygonNTetras[k]; ++j) {
1262 cellsConn[pointsPerCell[i] + j] = cellVertsIds[pointsPerCell[i] + j];
1266 TTK_PSORT(this->threadNumber_, cellVertsIds.begin(), cellVertsIds.end());
1267 const auto last = std::unique(cellVertsIds.begin(), cellVertsIds.end());
1268 cellVertsIds.erase(last, cellVertsIds.end());
1271 std::vector<SimplexId> vertId2PointsId(triangulation.getNumberOfCells());
1273 const auto noldpoints{npoints};
1274 npoints += cellVertsIds.size();
1275 ncells = noldcells + validTetraIds.size();
1278 outSeps2.pt.points_.resize(3 * npoints);
1279 auto points = &outSeps2.pt.points_[3 * noldpoints];
1280 outSeps2.cl.offsets_.resize(ncells + 1);
1281 outSeps2.cl.offsets_[0] = 0;
1282 auto cellsOff = &outSeps2.cl.offsets_[noldcells];
1283 outSeps2.cl.sourceIds_.resize(ncells);
1284 outSeps2.cl.separatrixIds_.resize(ncells);
1285 outSeps2.cl.separatrixTypes_.resize(ncells);
1286 outSeps2.cl.isOnBoundary_.resize(ncells);
1288#ifdef TTK_ENABLE_OPENMP
1289#pragma omp parallel for num_threads(threadNumber_)
1291 for(
size_t i = 0; i < cellVertsIds.size(); ++i) {
1293 triangulation.getTetraIncenter(cellVertsIds[i], &points[3 * i]);
1295 vertId2PointsId[cellVertsIds[i]] = i + noldpoints;
1298#ifdef TTK_ENABLE_OPENMP
1299#pragma omp parallel for num_threads(threadNumber_)
1301 for(
size_t i = 0; i < validTetraIds.size(); ++i) {
1302 const auto m = validTetraIds[i];
1303 const auto k = pointsPerCell[i];
1304 for(
SimplexId j = 0; j < polygonNTetras[m]; ++j) {
1305 cellsConn[k + j] = vertId2PointsId[cellsConn[k + j]];
1307 const auto l = i + noldcells;
1308 const auto n = polygonSepInfosIds[m];
1309 outSeps2.cl.sourceIds_[l] = sepSourceIds[n];
1310 outSeps2.cl.separatrixIds_[l] = sepIds[n];
1311 outSeps2.cl.separatrixTypes_[l] = 1;
1312 outSeps2.cl.isOnBoundary_[l] = sepOnBoundary[n];
1315 for(
size_t i = 0; i < validTetraIds.size(); ++i) {
1317 cellsOff[i + 1] = cellsOff[i] + polygonNTetras[validTetraIds[i]];
1320 outSeps2.pt.numberOfPoints_ = npoints;
1321 outSeps2.cl.numberOfCells_ = ncells;
1329 const std::vector<Separatrix> &separatrices,
1330 const std::vector<std::vector<SimplexId>> &separatricesSaddles,
1332 const triangulationType &triangulation)
const {
1350 const auto noldcells{ncells};
1355 std::vector<size_t> geomCellsBegId{ncells};
1358 for(
size_t i = 0; i < separatrices.size(); ++i) {
1359 const auto &sep = separatrices[i];
1360 ncells += sep.geometry_.size();
1361 geomCellsBegId.emplace_back(ncells);
1368 firstCellId + 3 * (ncells - noldcells));
1369 auto cellsOff = &outSeps2.
cl.
offsets_[noldcells];
1374 separatrixFunctionMaxima.resize(separatrixId + separatrices.size());
1375 separatrixFunctionMinima.resize(separatrixId + separatrices.size());
1379 std::vector<SimplexId> cellVertsIds(3 * (ncells - noldcells));
1381#ifdef TTK_ENABLE_OPENMP
1382#pragma omp parallel for num_threads(threadNumber_)
1384 for(
size_t i = 0; i < separatrices.size(); ++i) {
1385 const auto &sep = separatrices[i];
1386 const auto &sepGeom = sep.geometry_;
1387 const auto &sepSaddles = separatricesSaddles[i];
1388 const auto sepId = separatrixId + i;
1390 const char sepType = 2;
1393 const auto sepFuncMax
1396 if(!sepSaddles.empty()) {
1398 const auto minId = *std::min_element(
1399 sepSaddles.begin(), sepSaddles.end(),
1401 return offsets[discreteGradient_.getCellLowerVertex(
1402 Cell{1, a}, triangulation)]
1404 Cell{1, b}, triangulation)];
1411 const auto minId = *std::min_element(
1412 sepGeom.begin(), sepGeom.end(),
1413 [&triangulation, offsets,
this](
const Cell &a,
const Cell &b) {
1414 return offsets[discreteGradient_.getCellLowerVertex(a, triangulation)]
1415 < offsets[discreteGradient_.getCellLowerVertex(
1420 separatrixFunctionMaxima[sepId] = sepFuncMax;
1421 separatrixFunctionMinima[sepId] = sepFuncMin;
1424 const char onBoundary
1425 = (sepSaddles.empty()
1427 : std::count_if(sepSaddles.begin(), sepSaddles.end(),
1429 return triangulation.isEdgeOnBoundary(a);
1431 + triangulation.isTriangleOnBoundary(src.id_);
1433 for(
size_t j = 0; j < sepGeom.size(); ++j) {
1434 const auto &cell = sepGeom[j];
1438 triangulation.getTriangleVertex(cell.id_, 0, v0);
1439 triangulation.getTriangleVertex(cell.id_, 1, v1);
1440 triangulation.getTriangleVertex(cell.id_, 2, v2);
1443 const auto l = geomCellsBegId[i] + j;
1445 const auto m = l - noldcells;
1447 cellsConn[3 * m + 0] = v0;
1448 cellsConn[3 * m + 1] = v1;
1449 cellsConn[3 * m + 2] = v2;
1450 cellVertsIds[3 * m + 0] = v0;
1451 cellVertsIds[3 * m + 1] = v1;
1452 cellVertsIds[3 * m + 2] = v2;
1463 TTK_PSORT(this->threadNumber_, cellVertsIds.begin(), cellVertsIds.end());
1464 const auto last = std::unique(cellVertsIds.begin(), cellVertsIds.end());
1465 cellVertsIds.erase(last, cellVertsIds.end());
1468 std::vector<size_t> vertId2PointsId(triangulation.getNumberOfVertices());
1470 const auto noldpoints{npoints};
1471 npoints += cellVertsIds.size();
1472 outSeps2.pt.points_.resize(3 * npoints);
1473 auto points = &outSeps2.pt.points_[3 * noldpoints];
1475#ifdef TTK_ENABLE_OPENMP
1476#pragma omp parallel for num_threads(threadNumber_)
1478 for(
size_t i = 0; i < cellVertsIds.size(); ++i) {
1480 triangulation.getVertexPoint(
1481 cellVertsIds[i], points[3 * i + 0], points[3 * i + 1], points[3 * i + 2]);
1483 vertId2PointsId[cellVertsIds[i]] = i + noldpoints;
1486 const auto lastOffset = noldcells == 0 ? 0 : cellsOff[-1];
1488#ifdef TTK_ENABLE_OPENMP
1489#pragma omp parallel for num_threads(threadNumber_)
1491 for(
size_t i = 0; i < ncells - noldcells; ++i) {
1492 cellsOff[i] = 3 * i + lastOffset;
1493 cellsConn[3 * i + 0] = vertId2PointsId[cellsConn[3 * i + 0]];
1494 cellsConn[3 * i + 1] = vertId2PointsId[cellsConn[3 * i + 1]];
1495 cellsConn[3 * i + 2] = vertId2PointsId[cellsConn[3 * i + 2]];
1498 cellsOff[ncells - noldcells] = cellsOff[ncells - noldcells - 1] + 3;
1500 outSeps2.pt.numberOfPoints_ = npoints;
1501 outSeps2.cl.numberOfCells_ = ncells;
1512 const std::vector<SimplexId> &maxima,
1514 const triangulationType &triangulation)
const {
1516 if(morseSmaleManifold ==
nullptr) {
1517 this->
printErr(
"Could not compute ascending segmentation");
1523 const auto nVerts{triangulation.getNumberOfVertices()};
1524 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, -1);
1525 if(maxima.empty()) {
1530 const auto nCells{triangulation.getNumberOfCells()};
1531 std::vector<SimplexId> morseSmaleManifoldOnCells(nCells, -1);
1534 for(
const auto &
id : maxima) {
1536 morseSmaleManifoldOnCells[id] = nMax++;
1539 const auto dim{triangulation.getDimensionality()};
1542 auto getFaceStarNumber = &triangulationType::getTriangleStarNumber;
1543 auto getFaceStar = &triangulationType::getTriangleStar;
1546 getFaceStarNumber = &triangulationType::getEdgeStarNumber;
1547 getFaceStar = &triangulationType::getEdgeStar;
1548 }
else if(dim == 1) {
1550 getFaceStarNumber = &triangulationType::getVertexStarNumber;
1551 getFaceStar = &triangulationType::getVertexStar;
1555 std::vector<SimplexId> visited{};
1557 std::vector<uint8_t> isMarked(nCells, 0);
1559#ifdef TTK_ENABLE_OPENMP
1560#pragma omp parallel for num_threads(threadNumber_) firstprivate(visited)
1563 if(isMarked[i] == 1) {
1568 while(morseSmaleManifoldOnCells[curr] == -1) {
1569 if(isMarked[curr] == 1) {
1574 Cell{dim, curr}, triangulation,
true)};
1576 const auto nStars{(triangulation.*getFaceStarNumber)(paired)};
1578 (triangulation.*getFaceStar)(paired, j, next);
1584 visited.emplace_back(curr);
1591 for(
const auto el : visited) {
1592 morseSmaleManifoldOnCells[el] = morseSmaleManifoldOnCells[curr];
1597#ifdef TTK_ENABLE_OPENMP
1598#pragma omp parallel for num_threads(threadNumber_)
1601 if(triangulation.getVertexStarNumber(i) < 1) {
1606 triangulation.getVertexStar(i, 0, starId);
1608 morseSmaleManifold[i] = morseSmaleManifoldOnCells[starId];
1611 this->
printMsg(
" Ascending segmentation computed", 1.0, tm.getElapsedTime(),
1620 const std::vector<SimplexId> &minima,
1622 const triangulationType &triangulation)
const {
1624 if(morseSmaleManifold ==
nullptr) {
1625 this->
printErr(
"Could not compute descending segmentation");
1631 const auto nVerts{triangulation.getNumberOfVertices()};
1633 if(minima.size() == 1) {
1635 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, 0);
1639 std::fill(morseSmaleManifold, morseSmaleManifold + nVerts, -1);
1642 for(
const auto &cp : minima) {
1644 morseSmaleManifold[cp] = nMin++;
1647 std::vector<SimplexId> visited{};
1649#ifdef TTK_ENABLE_OPENMP
1650#pragma omp parallel for num_threads(threadNumber_) firstprivate(visited)
1653 if(morseSmaleManifold[i] != -1) {
1658 while(morseSmaleManifold[curr] == -1) {
1660 const auto pairedEdge{
1663 triangulation.getEdgeVertex(pairedEdge, 0, next);
1665 triangulation.getEdgeVertex(pairedEdge, 1, next);
1667 visited.emplace_back(curr);
1670 for(
const auto el : visited) {
1671 morseSmaleManifold[el] = morseSmaleManifold[curr];
1675 this->
printMsg(
" Descending segmentation computed", 1.0, tm.getElapsedTime(),
1685 const SimplexId *
const ascendingManifold,
1686 const SimplexId *
const descendingManifold,
1688 const triangulationType &triangulation)
const {
1690 if(ascendingManifold ==
nullptr || descendingManifold ==
nullptr
1691 || morseSmaleManifold ==
nullptr) {
1692 this->
printErr(
"Could not compute final segmentation");
1698 const size_t nVerts = triangulation.getNumberOfVertices();
1702#ifdef TTK_ENABLE_OPENMP
1703#pragma omp parallel for num_threads(threadNumber_)
1705 for(
size_t i = 0; i < nVerts; ++i) {
1706 const auto d = ascendingManifold[i];
1707 const auto a = descendingManifold[i];
1708 if(a == -1 || d == -1) {
1709 morseSmaleManifold[i] = -1;
1711 morseSmaleManifold[i] = a * numberOfMaxima + d;
1716 std::vector<SimplexId> sparseRegionIds(
1717 morseSmaleManifold, morseSmaleManifold + nVerts);
1721 this->
threadNumber_, sparseRegionIds.begin(), sparseRegionIds.end());
1722 const auto last = std::unique(sparseRegionIds.begin(), sparseRegionIds.end());
1723 sparseRegionIds.erase(last, sparseRegionIds.end());
1726 std::map<SimplexId, size_t> sparseToDenseRegionId{};
1728 for(
size_t i = 0; i < sparseRegionIds.size(); ++i) {
1729 sparseToDenseRegionId[sparseRegionIds[i]] = i;
1734#ifdef TTK_ENABLE_OPENMP
1735#pragma omp parallel for num_threads(threadNumber_)
1737 for(
size_t i = 0; i < nVerts; ++i) {
1738 morseSmaleManifold[i] = sparseToDenseRegionId[morseSmaleManifold[i]];
1741 this->
printMsg(
" Final segmentation computed", 1.0, tm.getElapsedTime(),
1749 const double persistenceThreshold,
1750 const dataType *
const scalars,
1752 const triangulationType &triangulation) {
1756 const auto dim{triangulation.getDimensionality()};
1758 this->
printWrn(
"Can't return saddle connectors without a 3D dataset");
1770 std::vector<PersPairType> dms_pairs{};
1771 dms.computePersistencePairs(dms_pairs, offsets, triangulation,
false,
true);
1776 const auto getPersistence
1777 = [
this, &triangulation, scalars](
const PersPairType &p) {
1779 Cell{2, p.death},
Cell{1, p.birth}, scalars, triangulation);
1783 const auto firstSadSadPair{std::distance(
1785 std::find_if(dms_pairs.begin(), dms_pairs.end(),
1786 [](
const auto &pair) { return pair.type == 1; }))};
1790 std::vector<bool> isVisited(triangulation.getNumberOfTriangles(),
false);
1791 std::vector<SimplexId> visitedTriangles{};
1796 std::vector<std::tuple<size_t, dataType>> pairs;
1797 for(
size_t i = firstSadSadPair; i < dms_pairs.size(); ++i) {
1798 const auto &pair{dms_pairs[i]};
1799 pairs.emplace_back(std::make_tuple(i, getPersistence(pair)));
1801 const auto comparePersistence
1802 = [](
const std::tuple<size_t, dataType> &pair1,
1803 const std::tuple<size_t, dataType> &pair2) {
1804 return std::get<1>(pair1) < std::get<1>(pair2);
1807 this->
threadNumber_, pairs.begin(), pairs.end(), comparePersistence);
1809 std::vector<std::tuple<dataType, SimplexId, SimplexId>> skippedPairsPers;
1812 for(
const auto &pairTup : pairs) {
1813 const auto &pairIndex = std::get<0>(pairTup);
1814 const auto &pair{dms_pairs[pairIndex]};
1815 const auto &pairPersistence = std::get<1>(pairTup);
1817 if(pair.type != 1 || pairPersistence > persistenceThreshold) {
1821 const Cell birth{1, pair.birth};
1822 const Cell death{2, pair.death};
1828 std::vector<Cell> vpath{};
1831 birth, death, isVisited, &vpath, triangulation, disableForkReversal);
1833 if(vpath.back() == death) {
1842 std::vector<bool> isVisitedUpdated(
1843 triangulation.getNumberOfTriangles(),
false);
1844 std::vector<SimplexId> visitedTrianglesUpdated{};
1845 VisitedMask maskUpdated{isVisitedUpdated, visitedTrianglesUpdated};
1848 birth, death, isVisitedUpdated,
nullptr, triangulation,
false,
true,
1853 vpath, triangulation,
true);
1854 this->
printMsg(
"Could not return saddle connector " + birth.to_string()
1855 +
" -> " + death.to_string()
1856 +
" without creating a cycle.",
1859 skippedPairsPers.emplace_back(
1860 std::make_tuple(pairPersistence, pair.birth, pair.death));
1865 this->
printMsg(
"Could not return saddle connector " + birth.to_string()
1866 +
" -> " + death.to_string(),
1869 skippedPairsPers.emplace_back(
1870 std::make_tuple(pairPersistence, pair.birth, pair.death));
1876 this->
threadNumber_, skippedPairsPers.begin(), skippedPairsPers.end());
1877 for(
unsigned int i = 0; i < skippedPairsPers.size(); ++i)
1878 this->
printMsg(std::to_string(i) +
" "
1879 + std::to_string(std::get<1>(skippedPairsPers[i])) +
" "
1880 + std::to_string(std::get<2>(skippedPairsPers[i])) +
" "
1881 + std::to_string(std::get<0>(skippedPairsPers[i])));
1884 this->
printMsg(
"Returned " + std::to_string(nReturned) +
" saddle connectors",
1885 1.0, tm.getElapsedTime(), this->threadNumber_);
1892 const triangulationType &triangulation,
1893 const double &persistenceThreshold,
1894 const bool &detectCriticalPoints) {
1896 std::vector<std::pair<SimplexId, char>> criticalPoints{};
1898 persistenceThreshold);
1900 if(detectCriticalPoints) {
1902 std::vector<Cell> criticalCells{};
1905 criticalPoints.resize(criticalCells.size());
1908 for(
size_t i = 0; i < criticalCells.size(); ++i) {
1909 const auto &c = criticalCells[i];
1919 const int numberOfDimensions
1921 std::vector<SimplexId> nDMTCriticalPoints(numberOfDimensions, 0);
1922 for(
const auto &c : criticalCells) {
1923 ++nDMTCriticalPoints[c.dim_];
1926 std::vector<SimplexId> nPLInteriorCriticalPoints(numberOfDimensions, 0);
1927 for(
const auto &cp : criticalPoints) {
1928 if(!triangulation.isVertexOnBoundary(cp.first)) {
1929 ++nPLInteriorCriticalPoints[cp.second];
1933 std::vector<std::vector<std::string>> rows(numberOfDimensions);
1934 for(
int i = 0; i < numberOfDimensions; ++i) {
1935 rows[i] = std::vector<std::string>{
1936 "#" + std::to_string(i) +
"-cell(s)",
1937 std::to_string(nDMTCriticalPoints[i]) +
" (with "
1938 + std::to_string(nPLInteriorCriticalPoints[i]) +
" interior PL)"};
1946 const bool allowBoundary =
true;
1948 std::vector<char> isPL;
1951 if(triangulation.getDimensionality() == 3) {
1953 allowBoundary, triangulation);
1957 "Gradient reversed", 1.0, t.
getElapsedTime(), this->threadNumber_);