TTK
Loading...
Searching...
No Matches
AbstractTriangulation.h
Go to the documentation of this file.
1
11
12#pragma once
13
14// base code includes
15#include <Cache.h>
16#include <Geometry.h>
17#include <Wrapper.h>
18
19#include <array>
20#include <ostream>
21#include <unordered_map>
22
23#ifdef TTK_ENABLE_KAMIKAZE
24#define TTK_TRIANGULATION_INTERNAL(NAME) NAME
25#else
26#define TTK_TRIANGULATION_INTERNAL(NAME) NAME##Internal
27#endif // TTK_ENABLE_KAMIKAZE
28
29#define ttkTemplateMacroCase(triangulationType, triangulationClass, call) \
30 case triangulationType: { \
31 using TTK_TT = triangulationClass; \
32 call; \
33 }; break
34
35#define ttkTemplateMacro(triangulationType, call) \
36 switch(triangulationType) { \
37 ttkTemplateMacroCase( \
38 ttk::Triangulation::Type::EXPLICIT, ttk::ExplicitTriangulation, call); \
39 ttkTemplateMacroCase( \
40 ttk::Triangulation::Type::IMPLICIT, ttk::ImplicitNoPreconditions, call); \
41 ttkTemplateMacroCase(ttk::Triangulation::Type::HYBRID_IMPLICIT, \
42 ttk::ImplicitWithPreconditions, call); \
43 ttkTemplateMacroCase( \
44 ttk::Triangulation::Type::PERIODIC, ttk::PeriodicNoPreconditions, call); \
45 ttkTemplateMacroCase(ttk::Triangulation::Type::HYBRID_PERIODIC, \
46 ttk::PeriodicWithPreconditions, call); \
47 ttkTemplateMacroCase( \
48 ttk::Triangulation::Type::COMPACT, ttk::CompactTriangulation, call); \
49 }
50
51namespace ttk {
52
53 // forward declaration of the ttk::dcg::DiscreteGradient class to
54 // give it access to the gradient cache (with a `friend`
55 // declaration)
56 namespace dcg {
57 class DiscreteGradient;
58 }
59
61
62 public:
64
66
71
74 int clear();
75
78 size_t footprint(size_t size = 0) const;
79
101 virtual inline int getCellEdge(const SimplexId &cellId,
102 const int &localEdgeId,
103 SimplexId &edgeId) const {
104
105#ifndef TTK_ENABLE_KAMIKAZE
106 // initialize output variable before early return
107 edgeId = -1;
108
110 return -1;
111#endif
112
113 if(getDimensionality() == 1)
114 return getCellNeighbor(cellId, localEdgeId, edgeId);
115
116 else if(getDimensionality() == 2)
117 return getTriangleEdgeInternal(cellId, localEdgeId, edgeId);
118
119 return getCellEdgeInternal(cellId, localEdgeId, edgeId);
120 }
121
138 virtual inline SimplexId getCellEdgeNumber(const SimplexId &cellId) const {
139#ifndef TTK_ENABLE_KAMIKAZE
141 return -1;
142#endif
143 if(getDimensionality() == 1)
144 return getCellNeighborNumber(cellId);
145
146 else if(getDimensionality() == 2)
147 return getTriangleEdgeNumber(cellId);
148
149 return getCellEdgeNumberInternal(cellId);
150 }
151
180 virtual inline const std::vector<std::vector<SimplexId>> *getCellEdges() {
181#ifndef TTK_ENABLE_KAMIKAZE
183 return nullptr;
184#endif
185 if(getDimensionality() == 1)
186 return getCellNeighbors();
187
188 else if(getDimensionality() == 2)
190
191 return getCellEdgesInternal();
192 }
193
213 virtual inline int getCellNeighbor(const SimplexId &cellId,
214 const int &localNeighborId,
215 SimplexId &neighborId) const {
216
217#ifndef TTK_ENABLE_KAMIKAZE
218 // initialize output variable before early return
219 neighborId = -1;
220
222 return -1;
223#endif
224 return getCellNeighborInternal(cellId, localNeighborId, neighborId);
225 }
226
241 virtual inline SimplexId
242 getCellNeighborNumber(const SimplexId &cellId) const {
243
244#ifndef TTK_ENABLE_KAMIKAZE
246 return -1;
247#endif
248 return getCellNeighborNumberInternal(cellId);
249 }
250
276 virtual inline const std::vector<std::vector<SimplexId>> *
278#ifndef TTK_ENABLE_KAMIKAZE
280 return nullptr;
281#endif
283 }
284
307 virtual inline int getCellTriangle(const SimplexId &cellId,
308 const int &localTriangleId,
309 SimplexId &triangleId) const {
310#ifndef TTK_ENABLE_KAMIKAZE
311 // initialize output variable before early return
312 triangleId = -1;
313
314 if(getDimensionality() == 1)
315 return -1;
316
318 return -2;
319#endif
320 if(getDimensionality() == 2)
321 return getCellNeighbor(cellId, localTriangleId, triangleId);
322
323 return getCellTriangleInternal(cellId, localTriangleId, triangleId);
324 }
325
343 virtual inline SimplexId
344 getCellTriangleNumber(const SimplexId &cellId) const {
345#ifndef TTK_ENABLE_KAMIKAZE
346 if(getDimensionality() == 1)
347 return -1;
348
350 return -2;
351#endif
352 if(getDimensionality() == 2)
353 return getCellNeighborNumber(cellId);
354
355 return getCellTriangleNumberInternal(cellId);
356 }
357
389 virtual inline const std::vector<std::vector<SimplexId>> *
391#ifndef TTK_ENABLE_KAMIKAZE
392 if(getDimensionality() == 1)
393 return nullptr;
394
396 return nullptr;
397#endif
398 if(getDimensionality() == 2)
399 return getCellNeighbors();
400
402 }
403
415 virtual inline int getCellVertex(const SimplexId &cellId,
416 const int &localVertexId,
417 SimplexId &vertexId) const {
418
419 return getCellVertexInternal(cellId, localVertexId, vertexId);
420 }
421
428 virtual inline SimplexId
429 getCellVertexNumber(const SimplexId &cellId) const {
430 return getCellVertexNumberInternal(cellId);
431 }
432
436 virtual inline int getDimensionality() const {
438 }
439
464 virtual inline const std::vector<std::array<SimplexId, 2>> *getEdges() {
465#ifndef TTK_ENABLE_KAMIKAZE
466 if(getDimensionality() == 1)
467 return nullptr;
468
470 return nullptr;
471#endif
472 return getEdgesInternal();
473 }
474
494 virtual inline int getEdgeLink(const SimplexId &edgeId,
495 const int &localLinkId,
496 SimplexId &linkId) const {
497#ifndef TTK_ENABLE_KAMIKAZE
498 // initialize output variable before early return
499 linkId = -1;
500
501 if(getDimensionality() == 1)
502 return -1;
503
505 return -2;
506#endif
507 return getEdgeLinkInternal(edgeId, localLinkId, linkId);
508 }
509
524 virtual inline SimplexId getEdgeLinkNumber(const SimplexId &edgeId) const {
525#ifndef TTK_ENABLE_KAMIKAZE
526 if(getDimensionality() == 1)
527 return -1;
528
530 return -2;
531#endif
532 return getEdgeLinkNumberInternal(edgeId);
533 }
534
558 virtual inline const std::vector<std::vector<SimplexId>> *getEdgeLinks() {
559#ifndef TTK_ENABLE_KAMIKAZE
560 if(getDimensionality() == 1)
561 return nullptr;
562
564 return nullptr;
565#endif
566 return getEdgeLinksInternal();
567 }
568
591 virtual inline int getEdgeStar(const SimplexId &edgeId,
592 const int &localStarId,
593 SimplexId &starId) const {
594
595#ifndef TTK_ENABLE_KAMIKAZE
596 // initialize output variable before early return
597 starId = -1;
598
599 if(getDimensionality() == 1)
600 return -1;
601
603 return -2;
604#endif
605 return getEdgeStarInternal(edgeId, localStarId, starId);
606 }
607
625 virtual inline SimplexId getEdgeStarNumber(const SimplexId &edgeId) const {
626#ifndef TTK_ENABLE_KAMIKAZE
627 if(getDimensionality() == 1)
628 return -1;
629
631 return -2;
632#endif
633 return getEdgeStarNumberInternal(edgeId);
634 }
635
664 virtual inline const std::vector<std::vector<SimplexId>> *getEdgeStars() {
665#ifndef TTK_ENABLE_KAMIKAZE
666 if(getDimensionality() == 1)
667 return nullptr;
668
670 return nullptr;
671#endif
672 return getEdgeStarsInternal();
673 }
674
694 virtual inline int getEdgeTriangle(const SimplexId &edgeId,
695 const int &localTriangleId,
696 SimplexId &triangleId) const {
697#ifndef TTK_ENABLE_KAMIKAZE
698 // initialize output variable before early return
699 triangleId = -1;
700
701 if(getDimensionality() == 1)
702 return -1;
703
705 return -2;
706#endif
707 if(getDimensionality() == 2)
708 return getEdgeStar(edgeId, localTriangleId, triangleId);
709
710 return getEdgeTriangleInternal(edgeId, localTriangleId, triangleId);
711 }
712
727 virtual inline SimplexId
728 getEdgeTriangleNumber(const SimplexId &edgeId) const {
729#ifndef TTK_ENABLE_KAMIKAZE
730 if(getDimensionality() == 1)
731 return -1;
732
734 return -2;
735#endif
736
737 if(getDimensionality() == 2)
738 return getEdgeStarNumber(edgeId);
739
740 return getEdgeTriangleNumberInternal(edgeId);
741 }
742
768 virtual inline const std::vector<std::vector<SimplexId>> *
770#ifndef TTK_ENABLE_KAMIKAZE
771 if(getDimensionality() == 1)
772 return nullptr;
773
775 return nullptr;
776#endif
777
778 if(getDimensionality() == 2)
779 return getEdgeStars();
780
782 }
783
801 virtual inline int getEdgeVertex(const SimplexId &edgeId,
802 const int &localVertexId,
803 SimplexId &vertexId) const {
804#ifndef TTK_ENABLE_KAMIKAZE
805 // initialize output variable before early return
806 vertexId = -1;
807
809 return -2;
810#endif
811 if(getDimensionality() == 1)
812 return getCellVertex(edgeId, localVertexId, vertexId);
813
814 return getEdgeVertexInternal(edgeId, localVertexId, vertexId);
815 }
816
819 inline int getEdgeVertexNumber(const SimplexId ttkNotUsed(edgeId)) const {
820 return 2;
821 }
822
830 virtual inline const std::array<SimplexId, 3> &getGridDimensions() const {
831 return this->gridDimensions_;
832 }
833
839 virtual inline SimplexId getNumberOfCells() const {
841 }
842
856 virtual inline SimplexId getNumberOfEdges() const {
857#ifndef TTK_ENABLE_KAMIKAZE
859 return -2;
860#endif
861 if(getDimensionality() == 1)
862 return getNumberOfCells();
863
865 }
866
880 virtual inline SimplexId getNumberOfTriangles() const {
881#ifndef TTK_ENABLE_KAMIKAZE
882 if(getDimensionality() == 1)
883 return -1;
884
886 return -2;
887#endif
888 if(getDimensionality() == 2)
889 return getNumberOfCells();
890
892 }
893
896 virtual inline SimplexId getNumberOfVertices() const {
898 }
899
924 virtual inline const std::vector<std::array<SimplexId, 3>> *getTriangles() {
925#ifndef TTK_ENABLE_KAMIKAZE
927 return nullptr;
928#endif
929 return getTrianglesInternal();
930 }
931
950 virtual inline int getTriangleEdge(const SimplexId &triangleId,
951 const int &localEdgeId,
952 SimplexId &edgeId) const {
953#ifndef TTK_ENABLE_KAMIKAZE
954 // initialize output variable before early return
955 edgeId = -1;
956
957 if(getDimensionality() == 1)
958 return -1;
959
961 return -2;
962#endif
963
964 return getTriangleEdgeInternal(triangleId, localEdgeId, edgeId);
965 }
966
981 virtual inline SimplexId
982 getTriangleEdgeNumber(const SimplexId &triangleId) const {
983#ifndef TTK_ENABLE_KAMIKAZE
984 if(getDimensionality() == 1)
985 return -1;
986
988 return -2;
989#endif
990
991 return getTriangleEdgeNumberInternal(triangleId);
992 }
993
1019 virtual inline const std::vector<std::vector<SimplexId>> *
1021#ifndef TTK_ENABLE_KAMIKAZE
1022 if(getDimensionality() == 1)
1023 return nullptr;
1024
1026 return nullptr;
1027#endif
1028
1029 return getTriangleEdgesInternal();
1030 }
1031
1051 virtual inline int getTriangleLink(const SimplexId &triangleId,
1052 const int &localLinkId,
1053 SimplexId &linkId) const {
1054#ifndef TTK_ENABLE_KAMIKAZE
1055 // initialize output variable before early return
1056 linkId = -1;
1057
1058 if(getDimensionality() != 3)
1059 return -1;
1060
1062 return -2;
1063#endif
1064 return getTriangleLinkInternal(triangleId, localLinkId, linkId);
1065 }
1066
1082 virtual inline SimplexId
1083 getTriangleLinkNumber(const SimplexId &triangleId) const {
1084#ifndef TTK_ENABLE_KAMIKAZE
1085 if(getDimensionality() != 3)
1086 return -1;
1087
1089 return -2;
1090#endif
1091 return getTriangleLinkNumberInternal(triangleId);
1092 }
1093
1118 virtual inline const std::vector<std::vector<SimplexId>> *
1120#ifndef TTK_ENABLE_KAMIKAZE
1121 if(getDimensionality() != 3)
1122 return nullptr;
1123
1125 return nullptr;
1126#endif
1127 return getTriangleLinksInternal();
1128 }
1129
1150 virtual inline int getTriangleStar(const SimplexId &triangleId,
1151 const int &localStarId,
1152 SimplexId &starId) const {
1153#ifndef TTK_ENABLE_KAMIKAZE
1154 // initialize output variable before early return
1155 starId = -1;
1156
1157 if(getDimensionality() != 3)
1158 return -1;
1159
1161 return -2;
1162#endif
1163 return getTriangleStarInternal(triangleId, localStarId, starId);
1164 }
1165
1180 virtual inline SimplexId
1181 getTriangleStarNumber(const SimplexId &triangleId) const {
1182#ifndef TTK_ENABLE_KAMIKAZE
1183 if(getDimensionality() != 3)
1184 return -1;
1185
1187 return -2;
1188#endif
1189 return getTriangleStarNumberInternal(triangleId);
1190 }
1191
1217 virtual inline const std::vector<std::vector<SimplexId>> *
1219#ifndef TTK_ENABLE_KAMIKAZE
1220 if(getDimensionality() != 3)
1221 return nullptr;
1222
1224 return nullptr;
1225#endif
1226 return getTriangleStarsInternal();
1227 }
1228
1246 virtual inline int getTriangleVertex(const SimplexId &triangleId,
1247 const int &localVertexId,
1248 SimplexId &vertexId) const {
1249#ifndef TTK_ENABLE_KAMIKAZE
1250 // initialize output variable before early return
1251 vertexId = -1;
1252
1253 if(getDimensionality() == 1)
1254 return -1;
1255
1257 return -2;
1258#endif
1259 if(getDimensionality() == 2)
1260 return getCellVertex(triangleId, localVertexId, vertexId);
1261
1262 return getTriangleVertexInternal(triangleId, localVertexId, vertexId);
1263 }
1264
1267 inline int
1269 return 3;
1270 }
1271
1292 virtual inline int getVertexEdge(const SimplexId &vertexId,
1293 const int &localEdgeId,
1294 SimplexId &edgeId) const {
1295
1296#ifndef TTK_ENABLE_KAMIKAZE
1297 // initialize output variable before early return
1298 edgeId = -1;
1299
1301 return -1;
1302#endif
1303 if(getDimensionality() == 1)
1304 return getVertexStar(vertexId, localEdgeId, edgeId);
1305
1306 return getVertexEdgeInternal(vertexId, localEdgeId, edgeId);
1307 }
1308
1323 virtual inline SimplexId
1324 getVertexEdgeNumber(const SimplexId &vertexId) const {
1325#ifndef TTK_ENABLE_KAMIKAZE
1327 return -1;
1328#endif
1329 if(getDimensionality() == 1)
1330 return getVertexStarNumber(vertexId);
1331
1332 return getVertexEdgeNumberInternal(vertexId);
1333 }
1334
1361 virtual inline const std::vector<std::vector<SimplexId>> *getVertexEdges() {
1362#ifndef TTK_ENABLE_KAMIKAZE
1364 return nullptr;
1365#endif
1366 if(getDimensionality() == 1)
1367 return getVertexStars();
1368
1369 return getVertexEdgesInternal();
1370 }
1371
1391 virtual inline int getVertexLink(const SimplexId &vertexId,
1392 const int &localLinkId,
1393 SimplexId &linkId) const {
1394#ifndef TTK_ENABLE_KAMIKAZE
1395 // initialize output variable before early return
1396 linkId = -1;
1397
1399 return -1;
1400#endif
1401 return getVertexLinkInternal(vertexId, localLinkId, linkId);
1402 }
1403
1418 virtual inline SimplexId
1419 getVertexLinkNumber(const SimplexId &vertexId) const {
1420#ifndef TTK_ENABLE_KAMIKAZE
1422 return -1;
1423#endif
1424 return getVertexLinkNumberInternal(vertexId);
1425 }
1426
1451 virtual inline const std::vector<std::vector<SimplexId>> *getVertexLinks() {
1452#ifndef TTK_ENABLE_KAMIKAZE
1454 return nullptr;
1455#endif
1456 return getVertexLinksInternal();
1457 }
1458
1475 virtual inline int getVertexNeighbor(const SimplexId &vertexId,
1476 const int &localNeighborId,
1477 SimplexId &neighborId) const {
1478#ifndef TTK_ENABLE_KAMIKAZE
1479 // initialize output variable before early return
1480 neighborId = -1;
1481
1483 return -1;
1484#endif
1485 return getVertexNeighborInternal(vertexId, localNeighborId, neighborId);
1486 }
1487
1499 virtual inline SimplexId
1500 getVertexNeighborNumber(const SimplexId &vertexId) const {
1501#ifndef TTK_ENABLE_KAMIKAZE
1503 return -1;
1504#endif
1505 return getVertexNeighborNumberInternal(vertexId);
1506 }
1507
1531 virtual inline const std::vector<std::vector<SimplexId>> *
1533#ifndef TTK_ENABLE_KAMIKAZE
1535 return nullptr;
1536#endif
1538 }
1539
1546 virtual inline int getVertexPoint(const SimplexId &vertexId,
1547 float &x,
1548 float &y,
1549 float &z) const {
1550 return getVertexPointInternal(vertexId, x, y, z);
1551 }
1552
1573 virtual inline int getVertexStar(const SimplexId &vertexId,
1574 const int &localStarId,
1575 SimplexId &starId) const {
1576#ifndef TTK_ENABLE_KAMIKAZE
1577 // initialize output variable before early return
1578 starId = -1;
1579
1581 return -1;
1582#endif
1583 return getVertexStarInternal(vertexId, localStarId, starId);
1584 }
1585
1600 virtual inline SimplexId
1601 getVertexStarNumber(const SimplexId &vertexId) const {
1602#ifndef TTK_ENABLE_KAMIKAZE
1604 return -1;
1605#endif
1606 return getVertexStarNumberInternal(vertexId);
1607 }
1608
1634 virtual inline const std::vector<std::vector<SimplexId>> *getVertexStars() {
1635#ifndef TTK_ENABLE_KAMIKAZE
1637 return nullptr;
1638#endif
1639 return getVertexStarsInternal();
1640 }
1641
1662 virtual inline int getVertexTriangle(const SimplexId &vertexId,
1663 const int &localTriangleId,
1664 SimplexId &triangleId) const {
1665#ifndef TTK_ENABLE_KAMIKAZE
1666 // initialize output variable before early return
1667 triangleId = -1;
1668
1669 if(getDimensionality() == 1)
1670 return -1;
1671
1673 return -2;
1674#endif
1675 if(getDimensionality() == 2)
1676 return getVertexStar(vertexId, localTriangleId, triangleId);
1677
1678 return getVertexTriangleInternal(vertexId, localTriangleId, triangleId);
1679 }
1680
1695 virtual inline SimplexId
1696 getVertexTriangleNumber(const SimplexId &vertexId) const {
1697#ifndef TTK_ENABLE_KAMIKAZE
1698 if(getDimensionality() == 1)
1699 return -1;
1700
1702 return -2;
1703#endif
1704
1705 if(getDimensionality() == 2)
1706 return getVertexStarNumber(vertexId);
1707
1708 return getVertexTriangleNumberInternal(vertexId);
1709 }
1710
1736 virtual inline const std::vector<std::vector<SimplexId>> *
1738#ifndef TTK_ENABLE_KAMIKAZE
1739 if(getDimensionality() == 1)
1740 return nullptr;
1741
1743 return nullptr;
1744#endif
1745 if(getDimensionality() == 2)
1746 return getVertexStars();
1747
1749 }
1750
1752 inline bool hasPeriodicBoundaries() const {
1754 }
1755
1775 virtual inline bool isEdgeOnBoundary(const SimplexId &edgeId) const {
1776#ifndef TTK_ENABLE_KAMIKAZE
1778 return false;
1779#endif
1780 return isEdgeOnBoundaryInternal(edgeId);
1781 }
1782
1785 virtual inline bool isEmpty() const {
1786 return true;
1787 }
1788
1792 virtual inline bool isManifold() const {
1793#ifndef TTK_ENABLE_KAMIKAZE
1794 if(!this->hasPreconditionedManifold())
1795 return false;
1796#endif
1797 return this->isManifold_;
1798 }
1799
1820 virtual inline int preconditionManifold() {
1821
1822 if(!this->hasPreconditionedManifold_) {
1824 this->hasPreconditionedManifold_ = true;
1825 }
1826 return 0;
1827 }
1828
1849 virtual inline bool
1850 isTriangleOnBoundary(const SimplexId &triangleId) const {
1851#ifndef TTK_ENABLE_KAMIKAZE
1853 return false;
1854#endif
1855 return isTriangleOnBoundaryInternal(triangleId);
1856 }
1857
1875 virtual inline bool isVertexOnBoundary(const SimplexId &vertexId) const {
1876#ifndef TTK_ENABLE_KAMIKAZE
1878 return false;
1879#endif
1880 return isVertexOnBoundaryInternal(vertexId);
1881 }
1882
1896 virtual inline int preconditionBoundaryEdges() {
1897
1902 }
1903 return 0;
1904 }
1905
1928
1942 virtual inline int preconditionBoundaryVertices() {
1943
1947 }
1948 return 0;
1949 }
1950
1966 virtual inline int preconditionCellEdges() {
1967
1970 if(getDimensionality() == 1)
1972
1975 }
1976 return 0;
1977 }
1978
1996 virtual inline int preconditionCellNeighbors() {
1997
2001 }
2002 return 0;
2003 }
2004
2022 virtual inline int preconditionCellTriangles() {
2023
2026
2027#ifndef TTK_ENABLE_KAMIKAZE
2028 if(getDimensionality() == 1)
2029 return -1;
2030#endif
2031 if(getDimensionality() == 2)
2033
2036 }
2037 return 0;
2038 }
2039
2057 virtual inline int preconditionEdges() {
2058
2062 }
2063 return 0;
2064 }
2065
2083 virtual inline int preconditionEdgeLinks() {
2084
2086#ifndef TTK_ENABLE_KAMIKAZE
2087 if(getDimensionality() == 1)
2088 return -1;
2089#endif
2093 }
2094 return 0;
2095 }
2096
2114 virtual inline int preconditionEdgeStars() {
2115
2117#ifndef TTK_ENABLE_KAMIKAZE
2118 if(getDimensionality() == 1)
2119 return -1;
2120#endif
2124 }
2125 return 0;
2126 }
2127
2145 virtual inline int preconditionEdgeTriangles() {
2146
2149
2150#ifndef TTK_ENABLE_KAMIKAZE
2151 if(getDimensionality() == 1)
2152 return -1;
2153#endif
2154
2155 if(getDimensionality() == 2)
2156 return preconditionEdgeStars();
2157
2161 }
2162
2163 return 0;
2164 }
2165
2183 virtual inline int preconditionTriangles() {
2184
2187
2188#ifndef TTK_ENABLE_KAMIKAZE
2189 if(getDimensionality() == 1)
2190 return -1;
2191#endif
2192 if(getDimensionality() == 2)
2193 return 0;
2194
2196 }
2197
2198 return 0;
2199 }
2200
2218 virtual inline int preconditionTriangleEdges() {
2219
2222
2223#ifndef TTK_ENABLE_KAMIKAZE
2224 if(getDimensionality() == 1)
2225 return -1;
2226#endif
2227 if(getDimensionality() == 2)
2228 return preconditionCellEdges();
2229
2233 }
2234 return 0;
2235 }
2236
2254 virtual inline int preconditionTriangleLinks() {
2255
2257#ifndef TTK_ENABLE_KAMIKAZE
2258 if(getDimensionality() != 3)
2259 return -2;
2260#endif
2264 }
2265 return 0;
2266 }
2267
2285 virtual inline int preconditionTriangleStars() {
2286
2288#ifndef TTK_ENABLE_KAMIKAZE
2289 if(getDimensionality() != 3)
2290 return -1;
2291#endif
2295 }
2296 return 0;
2297 }
2298
2316 virtual inline int preconditionVertexEdges() {
2317
2320
2321 if(getDimensionality() == 1)
2322 return preconditionVertexStars();
2323
2326 }
2327 return 0;
2328 }
2329
2347 virtual inline int preconditionVertexLinks() {
2348
2352 }
2353 return 0;
2354 }
2355
2373 virtual inline int preconditionVertexNeighbors() {
2374
2378 }
2379 return 0;
2380 }
2381
2399 virtual inline int preconditionVertexStars() {
2400
2404 }
2405 return 0;
2406 }
2407
2425 virtual inline int preconditionVertexTriangles() {
2426
2429
2430#ifndef TTK_ENABLE_KAMIKAZE
2431 if(getDimensionality() == 1)
2432 return -1;
2433#endif
2434
2435 if(getDimensionality() == 2)
2436 return preconditionVertexStars();
2437
2440 }
2441 return 0;
2442 }
2443
2447 inline int getEdgeIncenter(const SimplexId edgeId,
2448 float incenter[3]) const {
2449 std::array<SimplexId, 2> vertexId{};
2450 getEdgeVertex(edgeId, 0, vertexId[0]);
2451 getEdgeVertex(edgeId, 1, vertexId[1]);
2452
2453 std::array<float, 6> p{};
2454 getVertexPoint(vertexId[0], p[0], p[1], p[2]);
2455 getVertexPoint(vertexId[1], p[3], p[4], p[5]);
2456
2457 incenter[0] = 0.5 * p[0] + 0.5 * p[3];
2458 incenter[1] = 0.5 * p[1] + 0.5 * p[4];
2459 incenter[2] = 0.5 * p[2] + 0.5 * p[5];
2460
2461 return 0;
2462 }
2463
2467 inline int getTriangleIncenter(const SimplexId triangleId,
2468 float incenter[3]) const {
2469
2470 std::array<SimplexId, 3> vertexId{};
2471 if(getDimensionality() == 2) {
2472 getCellVertex(triangleId, 0, vertexId[0]);
2473 getCellVertex(triangleId, 1, vertexId[1]);
2474 getCellVertex(triangleId, 2, vertexId[2]);
2475 } else if(getDimensionality() == 3) {
2476 getTriangleVertex(triangleId, 0, vertexId[0]);
2477 getTriangleVertex(triangleId, 1, vertexId[1]);
2478 getTriangleVertex(triangleId, 2, vertexId[2]);
2479 }
2480
2481 std::array<float, 9> p{};
2482 getVertexPoint(vertexId[0], p[0], p[1], p[2]);
2483 getVertexPoint(vertexId[1], p[3], p[4], p[5]);
2484 getVertexPoint(vertexId[2], p[6], p[7], p[8]);
2485
2486 std::array<float, 3> d{};
2487 d[0] = Geometry::distance(&p[3], &p[6]);
2488 d[1] = Geometry::distance(&p[0], &p[6]);
2489 d[2] = Geometry::distance(&p[0], &p[3]);
2490 const float sum = d[0] + d[1] + d[2];
2491
2492 d[0] = d[0] / sum;
2493 d[1] = d[1] / sum;
2494 d[2] = d[2] / sum;
2495
2496 incenter[0] = d[0] * p[0] + d[1] * p[3] + d[2] * p[6];
2497 incenter[1] = d[0] * p[1] + d[1] * p[4] + d[2] * p[7];
2498 incenter[2] = d[0] * p[2] + d[1] * p[5] + d[2] * p[8];
2499
2500 return 0;
2501 }
2502
2507 inline int getTetraIncenter(const SimplexId tetraId,
2508 float incenter[3]) const {
2509 incenter[0] = 0.0f;
2510 incenter[1] = 0.0f;
2511 incenter[2] = 0.0f;
2512
2513 std::array<float, 3> p{};
2514 for(int i = 0; i < 4; ++i) {
2515 SimplexId triangleId{-1};
2516 getCellTriangle(tetraId, i, triangleId);
2517 getTriangleIncenter(triangleId, p.data());
2518 incenter[0] += p[0];
2519 incenter[1] += p[1];
2520 incenter[2] += p[2];
2521 }
2522
2523 incenter[0] /= 4.0f;
2524 incenter[1] /= 4.0f;
2525 incenter[2] /= 4.0f;
2526
2527 return 0;
2528 }
2529
2533 inline int getCellIncenter(const SimplexId cellid,
2534 const int dim,
2535 float incenter[3]) const {
2536 switch(dim) {
2537 case 0:
2538 getVertexPoint(cellid, incenter[0], incenter[1], incenter[2]);
2539 break;
2540 case 1:
2541 getEdgeIncenter(cellid, incenter);
2542 break;
2543 case 2:
2544 getTriangleIncenter(cellid, incenter);
2545 break;
2546 case 3:
2547 getTetraIncenter(cellid, incenter);
2548 break;
2549 }
2550 return 0;
2551 }
2552
2553 virtual inline int getCellVTKID(const int &ttkId, int &vtkId) const {
2554
2555#ifndef TTK_ENABLE_KAMIKAZE
2556 // initialize output variable before early return
2557 vtkId = -1;
2558#endif
2559 return getCellVTKIDInternal(ttkId, vtkId);
2560 }
2561
2562#ifdef TTK_ENABLE_MPI
2563
2564 // "vtkGhostType" on points & cells
2565
2566 inline void setVertexGhostArray(const unsigned char *const data) {
2567 this->vertexGhost_ = data;
2568 }
2569
2570 inline void setCellGhostArray(const unsigned char *const data) {
2571 this->cellGhost_ = data;
2572 }
2573
2584 virtual int preconditionGlobalBoundary() {
2585
2586 if(!hasPreconditionedGlobalBoundary_) {
2587 preconditionGlobalBoundaryInternal();
2588 hasPreconditionedGlobalBoundary_ = true;
2589 }
2590 return 0;
2591 }
2592
2593 virtual inline int preconditionGlobalBoundaryInternal() {
2594 return 0;
2595 }
2596
2597 // public preconditions methods (used in Triangulation/ttkAlgorithm)
2598 virtual int preconditionDistributedVertices() {
2599 return 0;
2600 }
2601 virtual int preconditionDistributedCells() {
2602 return 0;
2603 }
2604
2605 inline void setBoundingBox(const double *const bBox) {
2606 this->boundingBox_
2607 = {bBox[0], bBox[1], bBox[2], bBox[3], bBox[4], bBox[5]};
2608 }
2609
2610 // This method should be called to initialize and populate the
2611 // ghostCellsPerOwner_ and the remoteGhostCells_ attributes.
2612
2613 virtual int preconditionExchangeGhostCells() {
2614 return 0;
2615 }
2616
2617 // This method should be called to initialize and populate the
2618 // ghostVerticesPerOwner_ and the remoteGhostVertices_ attributes.
2619
2620 virtual int preconditionExchangeGhostVertices() {
2621 return 0;
2622 }
2623
2624 virtual int setVertexRankArray(const int *ttkNotUsed(rankArray)) {
2625 return 0;
2626 }
2627
2628 virtual int setCellRankArray(const int *ttkNotUsed(rankArray)) {
2629 return 0;
2630 }
2631
2632 virtual int preconditionEdgeRankArray() {
2633 return 0;
2634 }
2635
2636 virtual int preconditionTriangleRankArray() {
2637 return 0;
2638 }
2639
2640 // global <-> local id mappings
2641
2642 virtual inline SimplexId getVertexGlobalId(const SimplexId lvid) const {
2643 if(!ttk::isRunningWithMPI()) {
2644 return lvid;
2645 }
2646#ifndef TTK_ENABLE_KAMIKAZE
2647 if(this->getDimensionality() != 1 && this->getDimensionality() != 2
2648 && this->getDimensionality() != 3) {
2649 this->printErr("Only 1D, 2D and 3D datasets are supported");
2650 return -1;
2651 }
2652 if(!this->hasPreconditionedDistributedVertices_) {
2653 this->printErr("VertexGlobalId query without pre-process!");
2654 this->printErr(
2655 "Please call preconditionDistributedVertices() in a pre-process.");
2656 return -1;
2657 }
2658 if(lvid < 0 || lvid >= this->getNumberOfVertices()) {
2659 return -1;
2660 }
2661#endif // TTK_ENABLE_KAMIKAZE
2662 return this->getVertexGlobalIdInternal(lvid);
2663 }
2664 virtual inline SimplexId getVertexLocalId(const SimplexId gvid) const {
2665
2666 if(!ttk::isRunningWithMPI()) {
2667 return gvid;
2668 }
2669#ifndef TTK_ENABLE_KAMIKAZE
2670 if(this->getDimensionality() != 1 && this->getDimensionality() != 2
2671 && this->getDimensionality() != 3) {
2672 this->printErr("Only 1D, 2D and 3D datasets are supported");
2673 return -1;
2674 }
2675 if(!this->hasPreconditionedDistributedVertices_) {
2676 this->printErr("VertexLocalId query without pre-process!");
2677 this->printErr(
2678 "Please call preconditionDistributedVertices() in a pre-process.");
2679 return -1;
2680 }
2681#endif // TTK_ENABLE_KAMIKAZE
2682 return this->getVertexLocalIdInternal(gvid);
2683 }
2684
2685 virtual inline SimplexId getCellGlobalId(const SimplexId lcid) const {
2686#ifndef TTK_ENABLE_KAMIKAZE
2687 if(this->getDimensionality() != 1 && this->getDimensionality() != 2
2688 && this->getDimensionality() != 3) {
2689 this->printErr("Only 1D, 2D and 3D datasets are supported");
2690 return -1;
2691 }
2692 if(!this->hasPreconditionedDistributedCells_) {
2693 this->printErr("CellGlobalId query without pre-process!");
2694 this->printErr(
2695 "Please call preconditionDistributedCells() in a pre-process.");
2696 return -1;
2697 }
2698 if(lcid < 0 || lcid >= this->getNumberOfCells()) {
2699 return -1;
2700 }
2701#endif // TTK_ENABLE_KAMIKAZE
2702 if(!ttk::isRunningWithMPI()) {
2703 return lcid;
2704 }
2705 return this->getCellGlobalIdInternal(lcid);
2706 }
2707 virtual inline SimplexId getCellLocalId(const SimplexId gcid) const {
2708#ifndef TTK_ENABLE_KAMIKAZE
2709 if(this->getDimensionality() != 1 && this->getDimensionality() != 2
2710 && this->getDimensionality() != 3) {
2711 this->printErr("Only 1D, 2D and 3D datasets are supported");
2712 return -1;
2713 }
2714 if(!this->hasPreconditionedDistributedCells_) {
2715 this->printErr("CellLocalId query without pre-process!");
2716 this->printErr(
2717 "Please call preconditionDistributedCells() in a pre-process.");
2718 return -1;
2719 }
2720#endif // TTK_ENABLE_KAMIKAZE
2721 if(!ttk::isRunningWithMPI()) {
2722 return gcid;
2723 }
2724 return this->getCellLocalIdInternal(gcid);
2725 }
2726
2727 virtual inline SimplexId getEdgeGlobalId(const SimplexId leid) const {
2728 const auto dim{this->getDimensionality()};
2729#ifndef TTK_ENABLE_KAMIKAZE
2730 if(dim != 1 && dim != 2 && dim != 3) {
2731 this->printErr("Only 1D, 2D and 3D datasets are supported");
2732 return -1;
2733 }
2734 if(!this->hasPreconditionedDistributedEdges_) {
2735 this->printErr("EdgeGlobalId query without pre-process!");
2736 this->printErr(
2737 "Please call preconditionDistributedEdges() in a pre-process.");
2738 return -1;
2739 }
2740 if(leid < 0 || leid >= this->getNumberOfEdges()) {
2741 return -1;
2742 }
2743#endif // TTK_ENABLE_KAMIKAZE
2744 if(!ttk::isRunningWithMPI()) {
2745 return leid;
2746 }
2747 if(dim == 2 || dim == 3) {
2748 return this->getEdgeGlobalIdInternal(leid);
2749 } else if(dim == 1) {
2750 return this->getCellGlobalIdInternal(leid);
2751 }
2752 return -1;
2753 }
2754 virtual inline SimplexId getEdgeLocalId(const SimplexId geid) const {
2755 const auto dim{this->getDimensionality()};
2756#ifndef TTK_ENABLE_KAMIKAZE
2757 if(dim != 1 && dim != 2 && dim != 3) {
2758 this->printErr("Only 1D, 2D and 3D datasets are supported");
2759 return -1;
2760 }
2761 if(!this->hasPreconditionedDistributedEdges_) {
2762 this->printErr("EdgeLocalId query without pre-process!");
2763 this->printErr(
2764 "Please call preconditionDistributedEdges() in a pre-process.");
2765 return -1;
2766 }
2767#endif // TTK_ENABLE_KAMIKAZE
2768 if(!ttk::isRunningWithMPI()) {
2769 return geid;
2770 }
2771 if(dim == 2 || dim == 3) {
2772 return this->getEdgeLocalIdInternal(geid);
2773 } else if(dim == 1) {
2774 return this->getCellLocalIdInternal(geid);
2775 }
2776 return -1;
2777 }
2778
2779 virtual inline SimplexId getTriangleGlobalId(const SimplexId ltid) const {
2780 const auto dim{this->getDimensionality()};
2781#ifndef TTK_ENABLE_KAMIKAZE
2782 if(dim != 3 && dim != 2) {
2783 this->printErr("Only 2D and 3D datasets are supported");
2784 return -1;
2785 }
2786 if(!this->hasPreconditionedDistributedEdges_) {
2787 this->printErr("TriangleGlobalId query without pre-process!");
2788 this->printErr(
2789 "Please call preconditionDistributedTriangles() in a pre-process.");
2790 return -1;
2791 }
2792 if(ltid < 0 || ltid >= this->getNumberOfTriangles()) {
2793 return -1;
2794 }
2795#endif // TTK_ENABLE_KAMIKAZE
2796 if(!ttk::isRunningWithMPI()) {
2797 return ltid;
2798 }
2799 if(dim == 3) {
2800 return this->getTriangleGlobalIdInternal(ltid);
2801 } else if(dim == 2) {
2802 return this->getCellGlobalIdInternal(ltid);
2803 }
2804 return -1;
2805 }
2806 virtual inline SimplexId getTriangleLocalId(const SimplexId gtid) const {
2807 const auto dim{this->getDimensionality()};
2808#ifndef TTK_ENABLE_KAMIKAZE
2809 if(dim != 3 && dim != 2) {
2810 this->printErr("Only 2D and 3D datasets are supported");
2811 return -1;
2812 }
2813 if(!this->hasPreconditionedDistributedEdges_) {
2814 this->printErr("TriangleLocalId query without pre-process!");
2815 this->printErr(
2816 "Please call preconditionDistributedTriangles() in a pre-process.");
2817 return -1;
2818 }
2819#endif // TTK_ENABLE_KAMIKAZE
2820 if(!ttk::isRunningWithMPI()) {
2821 return gtid;
2822 }
2823 if(dim == 3) {
2824 return this->getTriangleLocalIdInternal(gtid);
2825 } else if(dim == 2) {
2826 return this->getCellLocalIdInternal(gtid);
2827 }
2828 return -1;
2829 }
2830
2831 virtual inline int getVertexRank(const SimplexId lvid) const {
2832
2833 if(!ttk::isRunningWithMPI()) {
2834 return 0;
2835 }
2836#ifndef TTK_ENABLE_KAMIKAZE
2837 if(this->getDimensionality() != 1 && this->getDimensionality() != 2
2838 && this->getDimensionality() != 3) {
2839 this->printErr("Only 1D, 2D and 3D datasets are supported");
2840 return -1;
2841 }
2842 if(!this->hasPreconditionedDistributedVertices_) {
2843 this->printErr("VertexRankId query without pre-process!");
2844 this->printErr(
2845 "Please call preconditionDistributedVertices() in a pre-process.");
2846 return -1;
2847 }
2848 if(lvid < 0 || lvid >= this->getNumberOfVertices()) {
2849 return -1;
2850 }
2851#endif // TTK_ENABLE_KAMIKAZE
2852 return this->getVertexRankInternal(lvid);
2853 }
2854
2855 virtual inline int getCellRank(const SimplexId lcid) const {
2856#ifndef TTK_ENABLE_KAMIKAZE
2857 if(this->getDimensionality() != 1 && this->getDimensionality() != 2
2858 && this->getDimensionality() != 3) {
2859 this->printErr("Only 1D, 2D and 3D datasets are supported");
2860 return -1;
2861 }
2862 if(!this->hasPreconditionedDistributedCells_) {
2863 this->printErr("CellRankId query without pre-process!");
2864 this->printErr(
2865 "Please call preconditionDistributedCells() in a pre-process.");
2866 return -1;
2867 }
2868 if(lcid < 0 || lcid >= this->getNumberOfCells()) {
2869 return -1;
2870 }
2871#endif // TTK_ENABLE_KAMIKAZE
2872 if(!ttk::isRunningWithMPI()) {
2873 return lcid;
2874 }
2875 return this->getCellRankInternal(lcid);
2876 }
2877
2878 virtual inline const std::vector<int> &getNeighborRanks() const {
2879 return this->neighborRanks_;
2880 }
2881 virtual inline std::vector<int> &getNeighborRanks() {
2882 return this->neighborRanks_;
2883 }
2884
2885 virtual inline std::map<int, int> &getNeighborsToId() {
2886 return this->neighborsToId_;
2887 }
2888
2889 virtual inline const std::map<int, int> &getNeighborsToId() const {
2890 return this->neighborsToId_;
2891 }
2892
2893 virtual inline const std::vector<std::array<ttk::SimplexId, 6>> &
2894 getNeighborVertexBBoxes() const {
2895 return this->neighborVertexBBoxes_;
2896 }
2897
2898 virtual inline const std::vector<std::vector<SimplexId>> &
2899 getGhostCellsPerOwner() const {
2900 if(!hasPreconditionedExchangeGhostCells_) {
2901 printErr("The ghostCellsOwner_ attribute has not been populated!");
2902 printErr(
2903 "Please call preconditionExchangeGhostCells in a pre-process.");
2904 }
2905 return this->ghostCellsPerOwner_;
2906 }
2907
2908 virtual inline const std::vector<std::vector<SimplexId>> &
2909 getRemoteGhostCells() const {
2910 if(!hasPreconditionedExchangeGhostCells_) {
2911 printErr("The remoteGhostCells_ attribute has not been populated!");
2912 printErr(
2913 "Please call preconditionExchangeGhostCells in a pre-process.");
2914 }
2915 return this->remoteGhostCells_;
2916 }
2917
2918 virtual inline const std::vector<std::vector<SimplexId>> &
2919 getGhostVerticesPerOwner() const {
2920 if(!hasPreconditionedExchangeGhostVertices_) {
2921 printErr(
2922 "The ghostVerticesPerOwner_ attribute has not been populated!");
2923 printErr(
2924 "Please call preconditionExchangeGhostVertices in a pre-process.");
2925 }
2926 return this->ghostVerticesPerOwner_;
2927 }
2928
2929 virtual inline const std::vector<std::vector<SimplexId>> &
2930 getRemoteGhostVertices() const {
2931 if(!hasPreconditionedExchangeGhostVertices_) {
2932 printErr("The remoteGhostVertices_ attribute has not been populated!");
2933 printErr(
2934 "Please call preconditionExchangeGhostVertices in a pre-process.");
2935 }
2936 return this->remoteGhostVertices_;
2937 }
2938
2939 virtual inline void setHasPreconditionedDistributedVertices(bool flag) {
2940 this->hasPreconditionedDistributedVertices_ = flag;
2941 }
2942
2943 virtual inline bool hasPreconditionedDistributedVertices() const {
2944 return this->hasPreconditionedDistributedVertices_;
2945 }
2946 virtual inline bool hasPreconditionedDistributedCells() const {
2947 return this->hasPreconditionedDistributedCells_;
2948 }
2949
2950 inline int getDistributedGlobalCellId(const SimplexId &localCellId,
2951 const int &cellDim,
2952 SimplexId &globalCellId) const {
2953 if(ttk::hasInitializedMPI()) {
2954 switch(cellDim) {
2955 case 0:
2956 globalCellId = this->getVertexGlobalIdInternal(localCellId);
2957 break;
2958 case 1:
2959 globalCellId = this->getEdgeGlobalIdInternal(localCellId);
2960 break;
2961 case 2:
2962 if(getDimensionality() == 2) {
2963 globalCellId = this->getCellGlobalIdInternal(localCellId);
2964 break;
2965 } else {
2966 globalCellId = this->getTriangleGlobalIdInternal(localCellId);
2967 break;
2968 }
2969 case 3: {
2970 globalCellId = this->getCellGlobalIdInternal(localCellId);
2971 break;
2972 }
2973 default:
2974 globalCellId = -1;
2975 break;
2976 }
2977 } else {
2978 globalCellId = localCellId;
2979 }
2980 return 0;
2981 }
2982
2983 inline bool isOrderArrayGlobal(const void *data) const {
2984#ifndef TTK_ENABLE_KAMIKAZE
2985 if(isOrderArrayGlobal_.find(data) != isOrderArrayGlobal_.end())
2986#endif
2987 return isOrderArrayGlobal_.at(data);
2988#ifndef TTK_ENABLE_KAMIKAZE
2989 else {
2990 printErr("No global array flag has been found for this scalar field");
2991 return false;
2992 }
2993#endif
2994 }
2995
2996 inline void setIsOrderArrayGlobal(const void *data, bool flag) {
2997 isOrderArrayGlobal_[data] = flag;
2998 }
2999
3000 protected:
3001 virtual inline SimplexId
3002 getVertexGlobalIdInternal(const SimplexId ttkNotUsed(lvid)) const {
3003 return -1;
3004 }
3005
3006 virtual inline SimplexId
3007 getVertexLocalIdInternal(const SimplexId ttkNotUsed(gvid)) const {
3008 return -1;
3009 }
3010
3011 // overriden in ImplicitTriangulation &
3012 // PeriodicImplicitTriangulation (where cellGid_ refers to
3013 // squares/cubes and not triangles/tetrahedron)
3014 virtual inline SimplexId
3015 getCellGlobalIdInternal(const SimplexId ttkNotUsed(lcid)) const {
3016 return -1;
3017 }
3018
3019 virtual inline SimplexId
3020 getCellLocalIdInternal(const SimplexId ttkNotUsed(gcid)) const {
3021 return -1;
3022 }
3023
3024 virtual inline SimplexId
3025 getEdgeGlobalIdInternal(const SimplexId ttkNotUsed(leid)) const {
3026 return -1;
3027 }
3028
3029 virtual inline SimplexId
3030 getEdgeLocalIdInternal(const SimplexId ttkNotUsed(geid)) const {
3031 return -1;
3032 }
3033
3034 virtual inline SimplexId
3035 getTriangleGlobalIdInternal(const SimplexId ttkNotUsed(ltid)) const {
3036 return -1;
3037 }
3038
3039 virtual inline SimplexId
3040 getTriangleLocalIdInternal(const SimplexId ttkNotUsed(gtid)) const {
3041 return -1;
3042 }
3043
3044 virtual inline int
3045 getVertexRankInternal(const SimplexId ttkNotUsed(lvid)) const {
3046 return -1;
3047 }
3048
3049 virtual inline int
3050 getCellRankInternal(const SimplexId ttkNotUsed(lcid)) const {
3051 return -1;
3052 }
3053
3054 // these protected methods should not be exposed since they are
3055 // called by their non-MPI-aware conterparts
3056 virtual inline bool
3057 isVertexOnGlobalBoundaryInternal(const SimplexId ttkNotUsed(lvid)) const {
3058 return false;
3059 }
3060 virtual inline bool
3061 isEdgeOnGlobalBoundaryInternal(const SimplexId ttkNotUsed(leid)) const {
3062 return false;
3063 }
3064 virtual inline bool isTriangleOnGlobalBoundaryInternal(
3065 const SimplexId ttkNotUsed(ltid)) const {
3066 return false;
3067 }
3068
3069#endif // TTK_ENABLE_MPI
3070
3071 protected:
3072 virtual int getCellEdgeInternal(const SimplexId &ttkNotUsed(cellId),
3073 const int &ttkNotUsed(localEdgeId),
3074 SimplexId &ttkNotUsed(edgeId)) const {
3075 return 0;
3076 }
3077
3078 virtual inline SimplexId
3080 return 0;
3081 }
3082
3083 virtual inline const std::vector<std::vector<SimplexId>> *
3085 return nullptr;
3086 }
3087
3088 virtual inline int
3090 const int &ttkNotUsed(localNeighborId),
3091 SimplexId &ttkNotUsed(neighborId)) const {
3092 return 0;
3093 }
3094
3095 virtual inline SimplexId
3097 return 0;
3098 }
3099
3100 virtual inline const std::vector<std::vector<SimplexId>> *
3102 return nullptr;
3103 }
3104
3105 virtual inline int
3107 const int &ttkNotUsed(localTriangleId),
3108 SimplexId &ttkNotUsed(triangleId)) const {
3109 return 0;
3110 }
3111
3112 virtual inline SimplexId
3114 return 0;
3115 }
3116
3117 virtual inline const std::vector<std::vector<SimplexId>> *
3119 return nullptr;
3120 }
3121
3122 virtual inline int
3124 const int &ttkNotUsed(localVertexId),
3125 SimplexId &ttkNotUsed(vertexId)) const {
3126 return 0;
3127 }
3128
3129 virtual inline SimplexId
3131 return 0;
3132 }
3133
3134 virtual inline int getDimensionalityInternal() const {
3135 return 0;
3136 }
3137
3138 virtual inline const std::vector<std::array<SimplexId, 2>> *
3140 return nullptr;
3141 }
3142
3143 virtual inline int
3145 const int &ttkNotUsed(localLinkId),
3146 SimplexId &ttkNotUsed(linkId)) const {
3147 return 0;
3148 }
3149
3150 virtual inline SimplexId
3152 return 0;
3153 }
3154
3155 virtual inline const std::vector<std::vector<SimplexId>> *
3157 return nullptr;
3158 }
3159
3160 virtual inline int
3162 const int &ttkNotUsed(localStarId),
3163 SimplexId &ttkNotUsed(starId)) const {
3164 return 0;
3165 }
3166
3167 virtual inline SimplexId
3169 return 0;
3170 }
3171
3172 virtual inline const std::vector<std::vector<SimplexId>> *
3174 return nullptr;
3175 }
3176
3177 virtual inline int
3179 const int &ttkNotUsed(localTriangleId),
3180 SimplexId &ttkNotUsed(triangleId)) const {
3181 return 0;
3182 }
3183
3184 virtual inline SimplexId
3186 return 0;
3187 }
3188
3189 virtual inline const std::vector<std::vector<SimplexId>> *
3191 return nullptr;
3192 }
3193
3194 virtual inline int
3196 const int &ttkNotUsed(localVertexId),
3197 SimplexId &ttkNotUsed(vertexId)) const {
3198 return 0;
3199 }
3200
3201 virtual inline SimplexId getNumberOfCellsInternal() const {
3202 return 0;
3203 }
3204
3205 virtual inline SimplexId getNumberOfEdgesInternal() const {
3206 return 0;
3207 }
3208
3210 return 0;
3211 }
3212
3214 return 0;
3215 }
3216
3217 virtual inline const std::vector<std::array<SimplexId, 3>> *
3219 return nullptr;
3220 }
3221
3222 virtual inline int
3224 const int &ttkNotUsed(localEdgeId),
3225 SimplexId &ttkNotUsed(edgeId)) const {
3226 return 0;
3227 }
3228
3230 const SimplexId &ttkNotUsed(triangleId)) const {
3231 return 0;
3232 }
3233
3234 virtual inline const std::vector<std::vector<SimplexId>> *
3236 return nullptr;
3237 }
3238
3239 virtual inline int
3241 const int &ttkNotUsed(localLinkId),
3242 SimplexId &ttkNotUsed(linkId)) const {
3243 return 0;
3244 }
3245
3247 const SimplexId &ttkNotUsed(triangleId)) const {
3248 return 0;
3249 }
3250
3251 virtual inline const std::vector<std::vector<SimplexId>> *
3253 return nullptr;
3254 }
3255
3256 virtual inline int
3258 const int &ttkNotUsed(localStarId),
3259 SimplexId &ttkNotUsed(starId)) const {
3260 return 0;
3261 }
3262
3264 const SimplexId &ttkNotUsed(triangleId)) const {
3265 return 0;
3266 }
3267
3268 virtual inline const std::vector<std::vector<SimplexId>> *
3270 return nullptr;
3271 }
3272
3273 virtual inline int
3275 const int &ttkNotUsed(localVertexId),
3276 SimplexId &ttkNotUsed(vertexId)) const {
3277 return 0;
3278 }
3279
3280 virtual inline int
3282 const int &ttkNotUsed(localEdgeId),
3283 SimplexId &ttkNotUsed(edgeId)) const {
3284 return 0;
3285 }
3286
3287 virtual inline SimplexId
3289 return 0;
3290 }
3291
3292 virtual inline const std::vector<std::vector<SimplexId>> *
3294 return nullptr;
3295 }
3296
3297 virtual inline int
3299 const int &ttkNotUsed(localLinkId),
3300 SimplexId &ttkNotUsed(linkId)) const {
3301 return 0;
3302 }
3303
3304 virtual inline SimplexId
3306 return 0;
3307 }
3308
3309 virtual inline const std::vector<std::vector<SimplexId>> *
3311 return nullptr;
3312 }
3313
3314 virtual inline int
3316 const int &ttkNotUsed(localNeighborId),
3317 SimplexId &ttkNotUsed(neighborId)) const {
3318 return 0;
3319 }
3320
3322 const SimplexId &ttkNotUsed(vertexId)) const {
3323 return 0;
3324 }
3325
3326 virtual inline const std::vector<std::vector<SimplexId>> *
3328 return nullptr;
3329 }
3330
3331 virtual inline int
3333 float &ttkNotUsed(x),
3334 float &ttkNotUsed(y),
3335 float &ttkNotUsed(z)) const {
3336 return 0;
3337 }
3338
3339 virtual inline int
3341 const int &ttkNotUsed(localStarId),
3342 SimplexId &ttkNotUsed(starId)) const {
3343 return 0;
3344 }
3345
3346 virtual inline SimplexId
3348 return 0;
3349 }
3350
3351 virtual inline const std::vector<std::vector<SimplexId>> *
3353 return nullptr;
3354 }
3355
3356 virtual inline int
3358 const int &ttkNotUsed(localTriangleId),
3359 SimplexId &ttkNotUsed(triangleId)) const {
3360 return 0;
3361 }
3362
3364 const SimplexId &ttkNotUsed(vertexId)) const {
3365 return 0;
3366 }
3367
3368 virtual inline const std::vector<std::vector<SimplexId>> *
3370 return nullptr;
3371 }
3372
3373 virtual inline bool
3375 return false;
3376 }
3377
3378 virtual inline bool isTriangleOnBoundaryInternal(
3379 const SimplexId &ttkNotUsed(triangleId)) const {
3380 return false;
3381 }
3382
3383 virtual inline bool
3385 return false;
3386 }
3387
3389
3390#ifndef TTK_ENABLE_KAMIKAZE
3392
3393 printMsg("BoundaryEdge query without pre-process!",
3395 printMsg("Please call preconditionBoundaryEdges() in a pre-process.",
3397 }
3398#endif
3399
3401 }
3402
3404#ifndef TTK_ENABLE_KAMIKAZE
3406
3407 printMsg("BoundaryTriangle query without pre-process!",
3409 printMsg(
3410 "Please call preconditionBoundaryTriangles() in a pre-process.",
3412 }
3413#endif
3415 }
3416
3418#ifndef TTK_ENABLE_KAMIKAZE
3420
3421 printMsg("BoundaryVertex query without pre-process!",
3423 printMsg("Please call preconditionBoundaryVertices() in a pre-process.",
3425 }
3426#endif
3428 }
3429
3430 inline bool hasPreconditionedCellEdges() const {
3431
3432#ifndef TTK_ENABLE_KAMIKAZE
3435
3436 printMsg("CellEdge query without pre-process!", debug::Priority::ERROR,
3437 debug::LineMode::NEW, std::cerr);
3438 printMsg("Please call preconditionCellEdges() in a pre-process.",
3440 return false;
3441 }
3442#endif
3443 return true;
3444 }
3445
3447#ifndef TTK_ENABLE_KAMIKAZE
3449
3450 printMsg("CellNeighbor query without pre-process!",
3452 printMsg("Please call preconditionCellNeighbors() in a pre-process.",
3454 }
3455#endif
3457 }
3458
3460#ifndef TTK_ENABLE_KAMIKAZE
3462 || ((getDimensionality() == 3)
3464
3465 printMsg("CellTriangle query without pre-process!",
3467 printMsg("Please call preconditionCellTriangles() in a pre-process.",
3469
3470 return false;
3471 }
3472#endif
3473 return true;
3474 }
3475
3476 inline bool hasPreconditionedEdgeLinks() const {
3477#ifndef TTK_ENABLE_KAMIKAZE
3479
3480 printMsg("EdgeLink query without pre-process!", debug::Priority::ERROR,
3481 debug::LineMode::NEW, std::cerr);
3482 printMsg("Please call preconditionEdgeLinks() in a pre-process.",
3484 }
3485#endif
3487 }
3488
3489 inline bool hasPreconditionedEdgeStars() const {
3490#ifndef TTK_ENABLE_KAMIKAZE
3492
3493 printMsg("EdgeStar query without pre-process!", debug::Priority::ERROR,
3494 debug::LineMode::NEW, std::cerr);
3495 printMsg("Please call preconditionEdgeStars() in a pre-process.",
3497 }
3498#endif
3500 }
3501
3503#ifndef TTK_ENABLE_KAMIKAZE
3505 || ((getDimensionality() == 3)
3507
3508 printMsg("EdgeTriangle query without pre-process!",
3510 printMsg("Please call preconditionEdgeTriangles() in a pre-process.",
3512
3513 return false;
3514 }
3515#endif
3516 return true;
3517 }
3518
3519 inline bool hasPreconditionedEdges() const {
3520#ifndef TTK_ENABLE_KAMIKAZE
3521 if((getDimensionality() != 1) && (!hasPreconditionedEdges_)) {
3522
3523 printMsg("Edge query without pre-process!", debug::Priority::ERROR,
3524 debug::LineMode::NEW, std::cerr);
3525 printMsg("Please call preconditionEdges() in a pre-process.",
3527
3528 return false;
3529 }
3530#endif
3531 return true;
3532 }
3533
3534 inline bool hasPreconditionedTriangles() const {
3535#ifndef TTK_ENABLE_KAMIKAZE
3537
3538 printMsg("Triangle query without pre-process!", debug::Priority::ERROR,
3539 debug::LineMode::NEW, std::cerr);
3540 printMsg("Please call preconditionTriangles() in a pre-process.",
3542
3543 return false;
3544 }
3545#endif
3546 return true;
3547 }
3548
3550#ifndef TTK_ENABLE_KAMIKAZE
3552 || ((getDimensionality() == 3)
3554
3555 printMsg("TriangleEdge query without pre-process!",
3557 printMsg("Please call preconditionTriangleEdges() in a pre-process.",
3559
3560 return false;
3561 }
3562#endif
3563 return true;
3564 }
3565
3567#ifndef TTK_ENABLE_KAMIKAZE
3569
3570 printMsg("TriangleLink query without pre-process!",
3572 printMsg("Please call preconditionTriangleLinks() in a pre-process.",
3574 }
3575#endif
3577 }
3578
3580#ifndef TTK_ENABLE_KAMIKAZE
3582
3583 printMsg("TriangleStar query without pre-process!",
3585 printMsg("Please call preconditionTriangleStars() in a pre-process.",
3587 }
3588#endif
3590 }
3591
3592 inline bool hasPreconditionedVertexEdges() const {
3593#ifndef TTK_ENABLE_KAMIKAZE
3596
3597 printMsg("VertexEdge query without pre-process!",
3599 printMsg("Please call preconditionVertexEdges() in a pre-process.",
3601
3602 return false;
3603 }
3604#endif
3605 return true;
3606 }
3607
3608 inline bool hasPreconditionedVertexLinks() const {
3609#ifndef TTK_ENABLE_KAMIKAZE
3611
3612 printMsg("VertexLink query without pre-process!",
3614 printMsg("Please call preconditionVertexLinks() in a pre-process.",
3616 }
3617#endif
3619 }
3620
3622#ifndef TTK_ENABLE_KAMIKAZE
3624
3625 printMsg("VertexNeighbor query without pre-process!",
3627 printMsg("Please call preconditionVertexNeighbors() in a pre-process.",
3629 }
3630#endif
3632 }
3633
3634 inline bool hasPreconditionedVertexStars() const {
3635#ifndef TTK_ENABLE_KAMIKAZE
3637
3638 printMsg("VertexStar query without pre-process!",
3640 printMsg("Please call preconditionVertexStars() in a pre-process.",
3642 }
3643#endif
3645 }
3646
3648#ifndef TTK_ENABLE_KAMIKAZE
3650 || ((getDimensionality() == 3)
3652
3653 printMsg("VertexTriangle query without pre-process!",
3655 printMsg("Please call preconditionVertexTriangles() in a pre-process.",
3657
3658 return false;
3659 }
3660#endif
3661 return true;
3662 }
3663
3664 inline bool hasPreconditionedManifold() const {
3665#ifndef TTK_ENABLE_KAMIKAZE
3666 if(!this->hasPreconditionedManifold_) {
3667 this->printMsg("isManifold query without pre-process!",
3669 this->printMsg("Please call preconditionManifold() in a pre-process.",
3671 }
3672#endif // TTK_ENABLE_KAMIKAZE
3673 return this->hasPreconditionedManifold_;
3674 }
3675
3676 // empty wrapping to VTK for now
3677 bool needsToAbort() override {
3678 return false;
3679 }
3680
3682 return 0;
3683 }
3684
3686 return 0;
3687 }
3688
3690 return 0;
3691 }
3692
3693 virtual inline int preconditionCellEdgesInternal() {
3694 return 0;
3695 }
3696
3698 return 0;
3699 }
3700
3702 return 0;
3703 }
3704
3705 virtual inline int preconditionEdgesInternal() {
3706 return 0;
3707 }
3708
3709 virtual inline int preconditionEdgeLinksInternal() {
3710 return 0;
3711 }
3712
3713 virtual inline int preconditionEdgeStarsInternal() {
3714 return 0;
3715 }
3716
3718 return 0;
3719 }
3720
3721 virtual inline int preconditionTrianglesInternal() {
3722 return 0;
3723 }
3724
3726 return 0;
3727 }
3728
3730 return 0;
3731 }
3732
3734 return 0;
3735 }
3736
3737 virtual inline int preconditionVertexEdgesInternal() {
3738 return 0;
3739 }
3740
3741 virtual inline int preconditionVertexLinksInternal() {
3742 return 0;
3743 }
3744
3746 return 0;
3747 }
3748
3749 virtual inline int preconditionVertexStarsInternal() {
3750 return 0;
3751 }
3752
3754 return 0;
3755 }
3756
3757 virtual inline int preconditionManifoldInternal() {
3758 return 0;
3759 }
3760
3761 virtual inline int getCellVTKIDInternal(const int &ttkId,
3762 int &vtkId) const {
3763#ifndef TTK_ENABLE_KAMIKAZE
3764 if(ttkId < 0) {
3765 return -1;
3766 }
3767#endif
3768 vtkId = ttkId;
3769 return 0;
3770 }
3771
3772 template <class itemType>
3773 size_t tableFootprint(const std::vector<itemType> &table,
3774 const std::string &tableName = "",
3775 std::ostream &stream = std::cout) const {
3776
3777 std::stringstream msg;
3778 if((table.size()) && (tableName.length())) {
3779 msg << tableName << ": " << table.size() * sizeof(itemType) << " bytes";
3780 printMsg(
3781 msg.str(), debug::Priority::INFO, debug::LineMode::NEW, stream);
3782 }
3783
3784 return table.size() * sizeof(itemType);
3785 }
3786
3787 template <class itemType>
3788 size_t tableTableFootprint(const std::vector<std::vector<itemType>> &table,
3789 const std::string &tableName = "",
3790 std::ostream &stream = std::cout) const;
3791
3792 int updateProgress(const float &ttkNotUsed(progress)) override {
3793 return 0;
3794 }
3795
3807
3808 bool isManifold_{true};
3809
3810 std::array<SimplexId, 3> gridDimensions_;
3811
3813 std::vector<std::array<SimplexId, 6>> tetraEdgeList_;
3814 std::vector<std::vector<SimplexId>> cellNeighborList_;
3815 std::vector<std::array<SimplexId, 4>> tetraTriangleList_;
3816 std::vector<std::vector<SimplexId>> edgeLinkList_;
3817 std::vector<std::array<SimplexId, 2>> edgeList_;
3818 std::vector<std::vector<SimplexId>> edgeStarList_;
3819 std::vector<std::vector<SimplexId>> edgeTriangleList_;
3820 std::vector<std::array<SimplexId, 3>> triangleList_;
3821 std::vector<std::array<SimplexId, 3>> triangleEdgeList_;
3822 std::vector<std::vector<SimplexId>> triangleLinkList_;
3823 std::vector<std::vector<SimplexId>> triangleStarList_;
3824 std::vector<std::vector<SimplexId>> vertexEdgeList_;
3825 std::vector<std::vector<SimplexId>> vertexLinkList_;
3826 std::vector<std::vector<SimplexId>> vertexNeighborList_;
3827 std::vector<std::vector<SimplexId>> vertexStarList_;
3828 std::vector<std::vector<SimplexId>> vertexTriangleList_;
3829
3830 // keep compatibility between getCellEdges(), getCellTriangles(),
3831 // getCellNeighbors() and getTriangleEdges()
3832 std::vector<std::vector<SimplexId>> cellEdgeVector_{};
3833 std::vector<std::vector<SimplexId>> cellTriangleVector_{};
3834 std::vector<std::vector<SimplexId>> triangleEdgeVector_{};
3835
3836#ifdef TTK_ENABLE_MPI
3837
3838 // precondition methods for distributed meshes
3839
3840 virtual int preconditionDistributedEdges() {
3841 return 0;
3842 }
3843 virtual int preconditionDistributedTriangles() {
3844 return 0;
3845 }
3846
3847 // "vtkGhostType" PointData array
3848 const unsigned char *vertexGhost_{};
3849 // "vtkGhostType" CellData array
3850 const unsigned char *cellGhost_{};
3851
3852 // list of neighboring ranks (sharing ghost cells to current rank)
3853 std::vector<int> neighborRanks_{};
3854 // list of ids in the neighborRanks_ vector of neighboring ranks
3855 std::map<int, int> neighborsToId_{};
3856 // global ids of (local) ghost cells per each MPI (neighboring) rank
3857 std::vector<std::vector<SimplexId>> ghostCellsPerOwner_{};
3858 // global ids of local (owned) cells that are ghost cells of other
3859 // (neighboring) ranks (per MPI rank)
3860 std::vector<std::vector<SimplexId>> remoteGhostCells_{};
3861 // global ids of (local) ghost vertices per each MPI (neighboring) rank
3862 std::vector<std::vector<SimplexId>> ghostVerticesPerOwner_{};
3863 // global ids of local (owned) vertices that are ghost cells of other
3864 // (neighboring) ranks (per MPI rank)
3865 std::vector<std::vector<SimplexId>> remoteGhostVertices_{};
3866 // hold the neighboring ranks vertex bounding boxes (metaGrid_ coordinates)
3867 std::vector<std::array<SimplexId, 6>> neighborVertexBBoxes_{};
3868 // hold the neighboring ranks cells bounding boxes (metaGrid_ coordinates)
3869 std::vector<std::array<double, 6>> neighborCellBBoxes_{};
3870 std::array<double, 6> boundingBox_{};
3871
3872 bool hasPreconditionedDistributedCells_{false};
3873 bool hasPreconditionedExchangeGhostCells_{false};
3874 bool hasPreconditionedDistributedEdges_{false};
3875 bool hasPreconditionedDistributedTriangles_{false};
3876 bool hasPreconditionedDistributedVertices_{false};
3877 bool hasPreconditionedExchangeGhostVertices_{false};
3878 bool hasPreconditionedGlobalBoundary_{false};
3879 std::map<const void *, bool> isOrderArrayGlobal_;
3880
3881#endif // TTK_ENABLE_MPI
3882
3883 // only ttk::dcg::DiscreteGradient should use what's defined below.
3885
3886#ifdef TTK_ENABLE_DCG_OPTIMIZE_MEMORY
3887 // might no be sufficient for Rips complexes, where a vertex can
3888 // have more than 128 neighbors
3889 using gradIdType = char;
3890#else
3892#endif
3893
3907 using gradientType = std::array<std::vector<gradIdType>, 6>;
3916 using gradientKeyType = std::pair<const void *, size_t>;
3917 /*
3918 * @brief Type for caching Discrete Gradient internal data structure.
3919 *
3920 * Uses the ttk::LRUCache with \ref gradientKeytype as key and
3921 * \ref gradientType as value types.
3922 */
3924 /*
3925 * @brief Access to the gradientCache_ mutable member variable
3926 *
3927 * \warning This should only be used by the
3928 * ttk::dcg::DiscreteGradient class.
3929 */
3931 return &this->gradientCache_;
3932 }
3933
3934 // store, for each triangulation object and per offset field, a
3935 // reference to the discrete gradient internal structure
3937 };
3938} // namespace ttk
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
AbstractTriangulation is an interface class that defines an interface for efficient traversal methods...
virtual const std::vector< std::vector< SimplexId > > * getCellNeighborsInternal()
virtual const std::vector< std::vector< SimplexId > > * getVertexNeighbors()
virtual int getVertexPoint(const SimplexId &vertexId, float &x, float &y, float &z) const
std::vector< std::vector< SimplexId > > vertexStarList_
virtual SimplexId getVertexEdgeNumberInternal(const SimplexId &ttkNotUsed(vertexId)) const
std::pair< const void *, size_t > gradientKeyType
Key type for gradientCacheType.
virtual int getEdgeLink(const SimplexId &edgeId, const int &localLinkId, SimplexId &linkId) const
virtual int getEdgeVertexInternal(const SimplexId &ttkNotUsed(edgeId), const int &ttkNotUsed(localVertexId), SimplexId &ttkNotUsed(vertexId)) const
virtual SimplexId getVertexTriangleNumber(const SimplexId &vertexId) const
virtual const std::vector< std::vector< SimplexId > > * getCellNeighbors()
virtual int getVertexTriangle(const SimplexId &vertexId, const int &localTriangleId, SimplexId &triangleId) const
virtual const std::vector< std::vector< SimplexId > > * getCellEdgesInternal()
virtual int getVertexEdge(const SimplexId &vertexId, const int &localEdgeId, SimplexId &edgeId) const
int getCellIncenter(const SimplexId cellid, const int dim, float incenter[3]) const
virtual SimplexId getEdgeTriangleNumber(const SimplexId &edgeId) const
virtual const std::array< SimplexId, 3 > & getGridDimensions() const
virtual bool isTriangleOnBoundary(const SimplexId &triangleId) const
virtual SimplexId getVertexTriangleNumberInternal(const SimplexId &ttkNotUsed(vertexId)) const
virtual const std::vector< std::vector< SimplexId > > * getVertexStars()
int getEdgeVertexNumber(const SimplexId ttkNotUsed(edgeId)) const
virtual const std::vector< std::vector< SimplexId > > * getEdgeStars()
virtual int getTriangleLinkInternal(const SimplexId &ttkNotUsed(triangleId), const int &ttkNotUsed(localLinkId), SimplexId &ttkNotUsed(linkId)) const
virtual SimplexId getNumberOfEdgesInternal() const
virtual const std::vector< std::vector< SimplexId > > * getVertexLinks()
virtual const std::vector< std::array< SimplexId, 2 > > * getEdgesInternal()
virtual const std::vector< std::vector< SimplexId > > * getVertexTriangles()
virtual int getTriangleStar(const SimplexId &triangleId, const int &localStarId, SimplexId &starId) const
virtual SimplexId getVertexNeighborNumber(const SimplexId &vertexId) const
bool hasPeriodicBoundaries() const
Returns true if the grid uses period boundary conditions.
virtual int getEdgeTriangle(const SimplexId &edgeId, const int &localTriangleId, SimplexId &triangleId) const
virtual const std::vector< std::vector< SimplexId > > * getTriangleEdgesInternal()
std::vector< std::array< SimplexId, 3 > > triangleList_
virtual int getCellEdgeInternal(const SimplexId &ttkNotUsed(cellId), const int &ttkNotUsed(localEdgeId), SimplexId &ttkNotUsed(edgeId)) const
virtual int getEdgeVertex(const SimplexId &edgeId, const int &localVertexId, SimplexId &vertexId) const
virtual const std::vector< std::vector< SimplexId > > * getVertexStarsInternal()
virtual int getEdgeLinkInternal(const SimplexId &ttkNotUsed(edgeId), const int &ttkNotUsed(localLinkId), SimplexId &ttkNotUsed(linkId)) const
AbstractTriangulation & operator=(AbstractTriangulation &&)=default
virtual SimplexId getNumberOfVerticesInternal() const
virtual const std::vector< std::array< SimplexId, 2 > > * getEdges()
virtual bool isVertexOnBoundaryInternal(const SimplexId &ttkNotUsed(vertexId)) const
virtual int getTriangleEdge(const SimplexId &triangleId, const int &localEdgeId, SimplexId &edgeId) const
virtual SimplexId getTriangleLinkNumberInternal(const SimplexId &ttkNotUsed(triangleId)) const
virtual SimplexId getNumberOfCellsInternal() const
std::vector< std::array< SimplexId, 4 > > tetraTriangleList_
int updateProgress(const float &ttkNotUsed(progress)) override
virtual SimplexId getTriangleLinkNumber(const SimplexId &triangleId) const
AbstractTriangulation & operator=(const AbstractTriangulation &)=default
std::vector< std::vector< SimplexId > > triangleLinkList_
virtual int getTriangleStarInternal(const SimplexId &ttkNotUsed(triangleId), const int &ttkNotUsed(localStarId), SimplexId &ttkNotUsed(starId)) const
virtual bool isVertexOnBoundary(const SimplexId &vertexId) const
virtual const std::vector< std::vector< SimplexId > > * getVertexEdgesInternal()
virtual const std::vector< std::vector< SimplexId > > * getTriangleLinks()
virtual SimplexId getCellEdgeNumber(const SimplexId &cellId) const
std::vector< std::vector< SimplexId > > vertexNeighborList_
std::vector< std::vector< SimplexId > > vertexLinkList_
virtual int getVertexTriangleInternal(const SimplexId &ttkNotUsed(vertexId), const int &ttkNotUsed(localTriangleId), SimplexId &ttkNotUsed(triangleId)) const
gradientCacheType * getGradientCacheHandler() const
virtual SimplexId getEdgeStarNumberInternal(const SimplexId &ttkNotUsed(edgeId)) const
virtual const std::vector< std::vector< SimplexId > > * getEdgeTrianglesInternal()
virtual bool isEdgeOnBoundary(const SimplexId &edgeId) const
virtual int getVertexNeighborInternal(const SimplexId &ttkNotUsed(vertexId), const int &ttkNotUsed(localNeighborId), SimplexId &ttkNotUsed(neighborId)) const
virtual int getCellNeighborInternal(const SimplexId &ttkNotUsed(cellId), const int &ttkNotUsed(localNeighborId), SimplexId &ttkNotUsed(neighborId)) const
AbstractTriangulation(const AbstractTriangulation &)=default
virtual const std::vector< std::vector< SimplexId > > * getEdgeTriangles()
virtual SimplexId getCellVertexNumberInternal(const SimplexId &ttkNotUsed(cellId)) const
virtual int getCellVTKIDInternal(const int &ttkId, int &vtkId) const
virtual const std::vector< std::array< SimplexId, 3 > > * getTriangles()
virtual int getTriangleLink(const SimplexId &triangleId, const int &localLinkId, SimplexId &linkId) const
virtual const std::vector< std::array< SimplexId, 3 > > * getTrianglesInternal()
virtual SimplexId getEdgeLinkNumberInternal(const SimplexId &ttkNotUsed(edgeId)) const
virtual int getVertexPointInternal(const SimplexId &ttkNotUsed(vertexId), float &ttkNotUsed(x), float &ttkNotUsed(y), float &ttkNotUsed(z)) const
virtual const std::vector< std::vector< SimplexId > > * getEdgeLinks()
virtual int getTriangleVertexInternal(const SimplexId &ttkNotUsed(triangleId), const int &ttkNotUsed(localVertexId), SimplexId &ttkNotUsed(vertexId)) const
std::vector< std::vector< SimplexId > > edgeLinkList_
virtual int getVertexLinkInternal(const SimplexId &ttkNotUsed(vertexId), const int &ttkNotUsed(localLinkId), SimplexId &ttkNotUsed(linkId)) const
std::vector< std::vector< SimplexId > > vertexEdgeList_
std::vector< std::vector< SimplexId > > triangleEdgeVector_
std::vector< std::vector< SimplexId > > triangleStarList_
virtual SimplexId getNumberOfCells() const
virtual SimplexId getCellTriangleNumber(const SimplexId &cellId) const
virtual SimplexId getNumberOfVertices() const
virtual SimplexId getVertexNeighborNumberInternal(const SimplexId &ttkNotUsed(vertexId)) const
virtual SimplexId getNumberOfTrianglesInternal() const
virtual int getEdgeStarInternal(const SimplexId &ttkNotUsed(edgeId), const int &ttkNotUsed(localStarId), SimplexId &ttkNotUsed(starId)) const
virtual SimplexId getCellNeighborNumber(const SimplexId &cellId) const
virtual int getDimensionality() const
virtual int getCellEdge(const SimplexId &cellId, const int &localEdgeId, SimplexId &edgeId) const
virtual int getEdgeStar(const SimplexId &edgeId, const int &localStarId, SimplexId &starId) const
virtual int getVertexNeighbor(const SimplexId &vertexId, const int &localNeighborId, SimplexId &neighborId) const
virtual SimplexId getCellEdgeNumberInternal(const SimplexId &ttkNotUsed(cellId)) const
std::vector< std::vector< SimplexId > > cellNeighborList_
virtual SimplexId getEdgeLinkNumber(const SimplexId &edgeId) const
virtual const std::vector< std::vector< SimplexId > > * getVertexEdges()
virtual const std::vector< std::vector< SimplexId > > * getCellTrianglesInternal()
virtual bool isTriangleOnBoundaryInternal(const SimplexId &ttkNotUsed(triangleId)) const
virtual int getCellTriangleInternal(const SimplexId &ttkNotUsed(cellId), const int &ttkNotUsed(localTriangleId), SimplexId &ttkNotUsed(triangleId)) const
virtual int getVertexStar(const SimplexId &vertexId, const int &localStarId, SimplexId &starId) const
virtual SimplexId getTriangleStarNumber(const SimplexId &triangleId) const
std::vector< std::array< SimplexId, 2 > > edgeList_
virtual int getCellVertex(const SimplexId &cellId, const int &localVertexId, SimplexId &vertexId) const
virtual SimplexId getCellVertexNumber(const SimplexId &cellId) const
int getTriangleIncenter(const SimplexId triangleId, float incenter[3]) const
virtual const std::vector< std::vector< SimplexId > > * getTriangleStars()
virtual const std::vector< std::vector< SimplexId > > * getTriangleLinksInternal()
virtual SimplexId getVertexLinkNumberInternal(const SimplexId &ttkNotUsed(vertexId)) const
std::vector< std::array< SimplexId, 6 > > tetraEdgeList_
std::array< SimplexId, 3 > gridDimensions_
virtual SimplexId getVertexStarNumberInternal(const SimplexId &ttkNotUsed(vertexId)) const
std::vector< std::vector< SimplexId > > edgeStarList_
virtual const std::vector< std::vector< SimplexId > > * getVertexTrianglesInternal()
virtual SimplexId getNumberOfTriangles() const
size_t tableFootprint(const std::vector< itemType > &table, const std::string &tableName="", std::ostream &stream=std::cout) const
virtual int getTriangleVertex(const SimplexId &triangleId, const int &localVertexId, SimplexId &vertexId) const
virtual int getTriangleEdgeInternal(const SimplexId &ttkNotUsed(triangleId), const int &ttkNotUsed(localEdgeId), SimplexId &ttkNotUsed(edgeId)) const
virtual SimplexId getVertexStarNumber(const SimplexId &vertexId) const
virtual int getCellVTKID(const int &ttkId, int &vtkId) const
virtual SimplexId getCellTriangleNumberInternal(const SimplexId &ttkNotUsed(cellId)) const
virtual SimplexId getCellNeighborNumberInternal(const SimplexId &ttkNotUsed(cellId)) const
virtual SimplexId getEdgeTriangleNumberInternal(const SimplexId &ttkNotUsed(edgeId)) const
virtual int preconditionBoundaryTrianglesInternal()
virtual bool isEdgeOnBoundaryInternal(const SimplexId &ttkNotUsed(edgeId)) const
int getTriangleVertexNumber(const SimplexId ttkNotUsed(triangleId)) const
virtual int getEdgeTriangleInternal(const SimplexId &ttkNotUsed(edgeId), const int &ttkNotUsed(localTriangleId), SimplexId &ttkNotUsed(triangleId)) const
int getEdgeIncenter(const SimplexId edgeId, float incenter[3]) const
virtual int getCellTriangle(const SimplexId &cellId, const int &localTriangleId, SimplexId &triangleId) const
virtual const std::vector< std::vector< SimplexId > > * getTriangleEdges()
virtual const std::vector< std::vector< SimplexId > > * getEdgeStarsInternal()
virtual int getVertexEdgeInternal(const SimplexId &ttkNotUsed(vertexId), const int &ttkNotUsed(localEdgeId), SimplexId &ttkNotUsed(edgeId)) const
virtual SimplexId getVertexLinkNumber(const SimplexId &vertexId) const
virtual const std::vector< std::vector< SimplexId > > * getVertexLinksInternal()
virtual SimplexId getVertexEdgeNumber(const SimplexId &vertexId) const
virtual SimplexId getTriangleEdgeNumberInternal(const SimplexId &ttkNotUsed(triangleId)) const
std::array< std::vector< gradIdType >, 6 > gradientType
Discrete gradient internal struct.
std::vector< std::vector< SimplexId > > vertexTriangleList_
virtual SimplexId getNumberOfEdges() const
std::vector< std::vector< SimplexId > > cellEdgeVector_
size_t footprint(size_t size=0) const
std::vector< std::vector< SimplexId > > cellTriangleVector_
size_t tableTableFootprint(const std::vector< std::vector< itemType > > &table, const std::string &tableName="", std::ostream &stream=std::cout) const
virtual const std::vector< std::vector< SimplexId > > * getCellTriangles()
int getTetraIncenter(const SimplexId tetraId, float incenter[3]) const
virtual const std::vector< std::vector< SimplexId > > * getTriangleStarsInternal()
virtual SimplexId getTriangleEdgeNumber(const SimplexId &triangleId) const
virtual int getCellNeighbor(const SimplexId &cellId, const int &localNeighborId, SimplexId &neighborId) const
virtual const std::vector< std::vector< SimplexId > > * getVertexNeighborsInternal()
virtual int getVertexLink(const SimplexId &vertexId, const int &localLinkId, SimplexId &linkId) const
std::vector< std::array< SimplexId, 3 > > triangleEdgeList_
virtual int getVertexStarInternal(const SimplexId &ttkNotUsed(vertexId), const int &ttkNotUsed(localStarId), SimplexId &ttkNotUsed(starId)) const
virtual SimplexId getEdgeStarNumber(const SimplexId &edgeId) const
virtual int getCellVertexInternal(const SimplexId &ttkNotUsed(cellId), const int &ttkNotUsed(localVertexId), SimplexId &ttkNotUsed(vertexId)) const
std::vector< std::vector< SimplexId > > edgeTriangleList_
virtual const std::vector< std::vector< SimplexId > > * getEdgeLinksInternal()
virtual SimplexId getTriangleStarNumberInternal(const SimplexId &ttkNotUsed(triangleId)) const
virtual const std::vector< std::vector< SimplexId > > * getCellEdges()
AbstractTriangulation(AbstractTriangulation &&)=default
virtual int getDimensionalityInternal() const
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
Wrapper class to wrap ttk code.
Definition Wrapper.h:14
TTK discreteGradient processing package.
T distance(const T *p0, const T *p1, const int &dimension=3)
Definition Geometry.cpp:362
The Topology ToolKit.
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)