36#include <boost/functional/hash.hpp>
37#include <boost/unordered_map.hpp>
38#include <boost/unordered_set.hpp>
47 std::vector<std::array<SimplexId, 2>> internalEdgeList_;
48 std::vector<std::array<SimplexId, 3>> internalTriangleList_;
49 boost::unordered_map<std::array<SimplexId, 2>,
SimplexId> internalEdgeMap_;
50 boost::unordered_map<std::array<SimplexId, 2>,
SimplexId> externalEdgeMap_;
51 boost::unordered_map<std::array<SimplexId, 3>,
SimplexId>
53 boost::unordered_map<std::array<SimplexId, 3>,
SimplexId>
56 std::vector<bool> boundaryVertices_;
57 std::vector<bool> boundaryEdges_;
58 std::vector<bool> boundaryTriangles_;
72 std::vector<std::array<SimplexId, 3>> triangleEdges_;
76 std::vector<std::array<SimplexId, 6>> tetraEdges_;
78 std::vector<std::array<SimplexId, 4>> tetraTriangles_;
92 enum class SIMPLEX_ID { EDGE_ID = 1, TRIANGLE_ID = 2 };
106 const void *pointSet,
107 const int *indexArray,
108 const bool &doublePrecision =
false) {
125#ifdef TTK_CELL_ARRAY_NEW
134 const auto &cellDimension = offset[1] - offset[0] - 1;
136 if(cellDimension < 0 || cellDimension > 3) {
137 this->
printErr(
"Unable to create triangulation for cells of "
138 "dimension 4 or higher ("
139 + std::to_string(cellDimension) +
").");
145#ifdef TTK_ENABLE_OPENMP
146#pragma omp parallel for num_threads(this->threadNumber_)
148 for(
SimplexId i = 0; i < cellNumber; i++) {
149 if(offset[i + 1] - offset[i] - 1 != cellDimension) {
150#ifdef TTK_ENABLE_OPENMP
151#pragma omp atomic write
158 this->
printErr(
"Unable to create triangulation for "
159 "inhomogeneous\ncell dimensions.");
171 reorderCells(vertexMap, cellNumber, connectivity, offset);
173 = std::make_shared<CellArray>(connectivity, offset, cellNumber);
203 cellArray, cellNumber, cellArray[0] - 1);
218 int reorderCells(
const std::vector<SimplexId> &vertexMap,
226 int reorderCells(
const std::vector<SimplexId> &vertexMap,
230 const int &localEdgeId,
233#ifndef TTK_ENABLE_KAMIKAZE
238 if(localEdgeId < 0) {
248 if(exnode->tetraEdges_.empty()) {
252 if(localEdgeId >= (
int)(exnode->tetraEdges_)[localCellId].size()) {
255 edgeId = (exnode->tetraEdges_)[localCellId][localEdgeId];
263#ifndef TTK_ENABLE_KAMIKAZE
272 inline const std::vector<std::vector<SimplexId>> *
278 if(exnode->tetraEdges_.empty()) {
281 for(
size_t i = 0; i < exnode->tetraEdges_.size(); i++) {
283 exnode->tetraEdges_.at(i).end());
292 const int &localNeighborId,
295#ifndef TTK_ENABLE_KAMIKAZE
300 if(localNeighborId < 0) {
310 if(exnode->cellNeighbors_.
empty()) {
314 if(localNeighborId >= exnode->cellNeighbors_.
size(localCellId)) {
317 neighborId = exnode->cellNeighbors_.
get(localCellId, localNeighborId);
323 const SimplexId &cellId)
const override {
325#ifndef TTK_ENABLE_KAMIKAZE
334 if(exnode->cellNeighbors_.
empty()) {
337 return exnode->cellNeighbors_.
size(localCellId);
340 inline const std::vector<std::vector<SimplexId>> *
344 std::vector<std::vector<SimplexId>> localCellNeighbors;
346 if(exnode->cellNeighbors_.
empty()) {
349 exnode->cellNeighbors_.
copyTo(localCellNeighbors);
351 localCellNeighbors.begin(),
352 localCellNeighbors.end());
358 const int &localTriangleId,
361#ifndef TTK_ENABLE_KAMIKAZE
366 if((localTriangleId < 0)
377 if(exnode->tetraTriangles_.empty()) {
380 triangleId = (exnode->tetraTriangles_)[localCellId][localTriangleId];
387#ifndef TTK_ENABLE_KAMIKAZE
396 inline const std::vector<std::vector<SimplexId>> *
402 if(exnode->tetraTriangles_.empty()) {
405 for(
size_t i = 0; i < exnode->tetraTriangles_.size(); i++) {
407 exnode->tetraTriangles_.at(i).begin(),
408 exnode->tetraTriangles_.at(i).end());
417 const int &localVertexId,
420#ifndef TTK_ENABLE_KAMIKAZE
425 if((localVertexId < 0)
426 || (localVertexId >=
cellArray_->getCellVertexNumber(cellId))) {
432 vertexId =
cellArray_->getCellVertex(cellId, localVertexId);
437 const SimplexId &cellId)
const override {
439#ifndef TTK_ENABLE_KAMIKAZE
444 return cellArray_->getCellVertexNumber(cellId);
451 inline const std::vector<std::array<SimplexId, 2>> *
456 if(exnode->internalEdgeList_.empty())
459 exnode->internalEdgeList_.end());
466 const int &localLinkId,
469#ifndef TTK_ENABLE_KAMIKAZE
474 if(localLinkId < 0) {
483 if(exnode->edgeLinks_.
empty()) {
487 if(localLinkId >= exnode->edgeLinks_.
size(localEdgeId)) {
490 linkId = exnode->edgeLinks_.
get(localEdgeId, localLinkId);
496 const SimplexId &edgeId)
const override {
498#ifndef TTK_ENABLE_KAMIKAZE
506 if(exnode->edgeLinks_.
empty()) {
509 return exnode->edgeLinks_.
size(localEdgeId);
512 inline const std::vector<std::vector<SimplexId>> *
516 std::vector<std::vector<SimplexId>> localEdgeLinks;
518 if(exnode->edgeLinks_.
empty()) {
521 exnode->edgeLinks_.
copyTo(localEdgeLinks);
523 edgeLinkList_.end(), localEdgeLinks.begin(), localEdgeLinks.end());
530 const int &localStarId,
533#ifndef TTK_ENABLE_KAMIKAZE
538 if(localStarId < 0) {
547 if(exnode->edgeStars_.
empty()) {
551 if(localStarId >= exnode->edgeStars_.
size(localEdgeId)) {
554 starId = exnode->edgeStars_.
get(localEdgeId, localStarId);
560 const SimplexId &edgeId)
const override {
562#ifndef TTK_ENABLE_KAMIKAZE
569 if(exnode->edgeStars_.
empty()) {
572 return exnode->edgeStars_.
size(localEdgeId);
575 inline const std::vector<std::vector<SimplexId>> *
579 std::vector<std::vector<SimplexId>> localEdgeStars;
581 if(exnode->edgeStars_.
empty()) {
584 exnode->edgeStars_.
copyTo(localEdgeStars);
586 edgeStarList_.end(), localEdgeStars.begin(), localEdgeStars.end());
592 const int &localTriangleId,
595#ifndef TTK_ENABLE_KAMIKAZE
600 if(localTriangleId < 0) {
609 if(exnode->edgeTriangles_.
empty()) {
613 if(localTriangleId >= exnode->edgeTriangles_.
size(localEdgeId)) {
616 triangleId = exnode->edgeTriangles_.
get(localEdgeId, localTriangleId);
624#ifndef TTK_ENABLE_KAMIKAZE
632 if(exnode->edgeTriangles_.
empty()) {
635 return exnode->edgeTriangles_.
size(localEdgeId);
638 inline const std::vector<std::vector<SimplexId>> *
642 std::vector<std::vector<SimplexId>> localEdgeTriangles;
644 if(exnode->edgeTriangles_.
empty()) {
647 exnode->edgeTriangles_.
copyTo(localEdgeTriangles);
649 localEdgeTriangles.begin(),
650 localEdgeTriangles.end());
656 const int &localVertexId,
659#ifndef TTK_ENABLE_KAMIKAZE
664 if((localVertexId != 0) && (localVertexId != 1)) {
673 if(exnode->internalEdgeList_.empty()) {
678 vertexId = exnode->internalEdgeList_.at(localEdgeId)[1];
680 vertexId = exnode->internalEdgeList_.at(localEdgeId)[0];
692#ifndef TTK_ENABLE_KAMIKAZE
702#ifndef TTK_ENABLE_KAMIKAZE
715 inline const std::vector<std::array<SimplexId, 3>> *
729 if(exnode->internalTriangleList_.empty()) {
733 exnode->internalTriangleList_.begin(),
734 exnode->internalTriangleList_.end());
741 const int &localEdgeId,
744#ifndef TTK_ENABLE_KAMIKAZE
749 if((localEdgeId < 0) || (localEdgeId > 2)) {
759 if(exnode->triangleEdges_.empty()) {
762 edgeId = (exnode->triangleEdges_)[localTriangleId][localEdgeId];
767 const SimplexId &triangleId)
const override {
768#ifndef TTK_ENABLE_KAMIKAZE
777 inline const std::vector<std::vector<SimplexId>> *
783 if(exnode->triangleEdges_.empty()) {
786 for(
size_t i = 0; i < exnode->triangleEdges_.size(); i++) {
788 exnode->triangleEdges_.at(i).begin(),
789 exnode->triangleEdges_.at(i).end());
798 const int &localLinkId,
801#ifndef TTK_ENABLE_KAMIKAZE
806 if(localLinkId < 0) {
816 if(exnode->triangleLinks_.
empty()) {
820 if(localLinkId >= exnode->triangleLinks_.
size(localTriangleId)) {
823 linkId = exnode->triangleLinks_.
get(localTriangleId, localLinkId);
829 const SimplexId &triangleId)
const override {
831#ifndef TTK_ENABLE_KAMIKAZE
840 if(exnode->triangleLinks_.
empty()) {
843 return exnode->triangleLinks_.
size(localTriangleId);
846 inline const std::vector<std::vector<SimplexId>> *
850 std::vector<std::vector<SimplexId>> localTriangleLinks;
852 if(exnode->triangleLinks_.
empty()) {
855 exnode->triangleLinks_.
copyTo(localTriangleLinks);
857 localTriangleLinks.begin(),
858 localTriangleLinks.end());
865 const int &localStarId,
868#ifndef TTK_ENABLE_KAMIKAZE
874 if(localStarId < 0) {
884 if(exnode->triangleStars_.
empty()) {
888 if(localStarId >= exnode->triangleStars_.
size(localTriangleId)) {
891 starId = exnode->triangleStars_.
get(localTriangleId, localStarId);
897 const SimplexId &triangleId)
const override {
899#ifndef TTK_ENABLE_KAMIKAZE
909 if(exnode->triangleStars_.
empty()) {
912 return exnode->triangleStars_.
size(localTriangleId);
915 inline const std::vector<std::vector<SimplexId>> *
919 std::vector<std::vector<SimplexId>> localTriangleStars;
921 if(exnode->triangleStars_.
empty()) {
925 localTriangleStars.begin(),
926 localTriangleStars.end());
932 const int &localVertexId,
935#ifndef TTK_ENABLE_KAMIKAZE
940 if((localVertexId < 0) || (localVertexId > 2)) {
950 if(exnode->internalTriangleList_.empty()) {
954 = exnode->internalTriangleList_.at(localTriangleId)[localVertexId];
959 const int &localEdgeId,
962#ifndef TTK_ENABLE_KAMIKAZE
967 if(localEdgeId < 0) {
976 if(exnode->vertexEdges_.
empty()) {
979 if(localEdgeId >= exnode->vertexEdges_.
size(localVertexId)) {
982 edgeId = exnode->vertexEdges_.
get(localVertexId, localEdgeId);
990#ifndef TTK_ENABLE_KAMIKAZE
998 if(exnode->vertexEdges_.
empty()) {
1001 return exnode->vertexEdges_.
size(localVertexId);
1004 inline const std::vector<std::vector<SimplexId>> *
1008 std::vector<std::vector<SimplexId>> localVertexEdges;
1010 if(exnode->vertexEdges_.
empty()) {
1013 exnode->vertexEdges_.
copyTo(localVertexEdges);
1015 localVertexEdges.end());
1023 const int &localLinkId,
1026#ifndef TTK_ENABLE_KAMIKAZE
1031 if(localLinkId < 0) {
1040 if(exnode->vertexLinks_.
empty()) {
1043 if(localLinkId >= exnode->vertexLinks_.
size(localVertexId)) {
1046 linkId = exnode->vertexLinks_.
get(localVertexId, localLinkId);
1052 const SimplexId &vertexId)
const override {
1054#ifndef TTK_ENABLE_KAMIKAZE
1062 if(exnode->vertexLinks_.
empty()) {
1065 return exnode->vertexLinks_.
size(localVertexId);
1068 inline const std::vector<std::vector<SimplexId>> *
1072 std::vector<std::vector<SimplexId>> localVertexLinks;
1074 if(exnode->vertexLinks_.
empty()) {
1077 exnode->vertexLinks_.
copyTo(localVertexLinks);
1079 localVertexLinks.end());
1086 const int &localNeighborId,
1088#ifndef TTK_ENABLE_KAMIKAZE
1093 if(localNeighborId < 0) {
1102 if(exnode ==
nullptr) {
1105 if(exnode->vertexNeighbors_.
empty()) {
1108 if(localNeighborId >= exnode->vertexNeighbors_.
size(localVertexId)) {
1112 = exnode->vertexNeighbors_.
get(localVertexId, localNeighborId);
1118 const SimplexId &vertexId)
const override {
1120#ifndef TTK_ENABLE_KAMIKAZE
1128 if(exnode->vertexNeighbors_.
empty()) {
1131 return exnode->vertexNeighbors_.
size(localVertexId);
1134 inline const std::vector<std::vector<SimplexId>> *
1138 std::vector<std::vector<SimplexId>> localVertexNeighbors;
1140 if(exnode->vertexNeighbors_.
empty()) {
1143 exnode->vertexNeighbors_.
copyTo(localVertexNeighbors);
1145 localVertexNeighbors.begin(),
1146 localVertexNeighbors.end());
1152 const SimplexId &vertexId,
float &x,
float &y,
float &z)
const override {
1154#ifndef TTK_ENABLE_KAMIKAZE
1160 x = ((
const double *)
pointSet_)[3 * vertexId];
1161 y = ((
const double *)
pointSet_)[3 * vertexId + 1];
1162 z = ((
const double *)
pointSet_)[3 * vertexId + 2];
1164 x = ((
const float *)
pointSet_)[3 * vertexId];
1165 y = ((
const float *)
pointSet_)[3 * vertexId + 1];
1166 z = ((
const float *)
pointSet_)[3 * vertexId + 2];
1174 const int &localStarId,
1177#ifndef TTK_ENABLE_KAMIKAZE
1182 if(localStarId < 0) {
1191 if(exnode->vertexStars_.
empty()) {
1194 if(localStarId >= exnode->vertexStars_.
size(localVertexId)) {
1197 starId = exnode->vertexStars_.
get(localVertexId, localStarId);
1203 const SimplexId &vertexId)
const override {
1205#ifndef TTK_ENABLE_KAMIKAZE
1213 if(exnode->vertexStars_.
empty()) {
1216 return exnode->vertexStars_.
size(localVertexId);
1219 inline const std::vector<std::vector<SimplexId>> *
1223 std::vector<std::vector<SimplexId>> localVertexStars;
1225 if(exnode->vertexStars_.
empty()) {
1228 exnode->vertexStars_.
copyTo(localVertexStars);
1230 localVertexStars.end());
1236 const int &localTriangleId,
1239#ifndef TTK_ENABLE_KAMIKAZE
1244 if(localTriangleId < 0) {
1253 if(exnode->vertexTriangles_.
empty()) {
1256 if(localTriangleId >= exnode->vertexTriangles_.
size(localVertexId)) {
1260 = exnode->vertexTriangles_.
get(localVertexId, localTriangleId);
1266 const SimplexId &vertexId)
const override {
1268#ifndef TTK_ENABLE_KAMIKAZE
1276 if(exnode->vertexTriangles_.
empty()) {
1279 return exnode->vertexTriangles_.
size(localVertexId);
1282 inline const std::vector<std::vector<SimplexId>> *
1286 std::vector<std::vector<SimplexId>> localVertexTriangles;
1288 if(exnode->vertexTriangles_.
empty()) {
1291 exnode->vertexTriangles_.
copyTo(localVertexTriangles);
1293 localVertexTriangles.begin(),
1294 localVertexTriangles.end());
1300 const SimplexId &edgeId)
const override {
1301#ifndef TTK_ENABLE_KAMIKAZE
1309 return (exnode->boundaryEdges_)[localedgeId];
1317 const SimplexId &triangleId)
const override {
1321#ifndef TTK_ENABLE_KAMIKAZE
1330 return (exnode->boundaryTriangles_)[localtriangleId];
1334 const SimplexId &vertexId)
const override {
1335#ifndef TTK_ENABLE_KAMIKAZE
1343 return (exnode->boundaryVertices_)[localVertexId];
1377#ifndef TTK_ENABLE_KAMIKAZE
1392 std::vector<SimplexId> edgeCount(
nodeNumber_ + 1);
1394#ifdef TTK_ENABLE_OPENMP
1395#pragma omp parallel for num_threads(threadNumber_)
1405 this->
printMsg(
"Edges preconditioned in "
1425#ifndef TTK_ENABLE_KAMIKAZE
1441 std::vector<SimplexId> triangleCount(
nodeNumber_ + 1);
1442#ifdef TTK_ENABLE_OPENMP
1443#pragma omp parallel for num_threads(threadNumber_)
1454 this->
printMsg(
"Triangles preconditioned in "
1483 this->
printErr(
"Unsupported dimension for vertex link precondition");
1495#ifndef TTK_ENABLE_KAMIKAZE
1565 const std::vector<SimplexId> *intervals =
nullptr;
1567 if(
idType == SIMPLEX_ID::EDGE_ID) {
1569 }
else if(
idType == SIMPLEX_ID::TRIANGLE_ID) {
1575 const std::vector<SimplexId>::const_iterator low
1576 = lower_bound(intervals->begin(), intervals->end(),
id);
1577 return (low - intervals->begin());
1586#ifdef TTK_ENABLE_OPENMP
1587 threadId = omp_get_thread_num();
1593 if(
caches_[threadId].back().nid == reservedId) {
1600 caches_[threadId].emplace_front(localCluster);
1610 bool computeInternalEdgeList,
1611 bool computeInternalEdgeMap)
const;
1621 bool computeInternalTriangleList,
1622 bool computeInternalTriangleMap)
const;
1734 mutable std::vector<std::list<ImplicitCluster>>
caches_;
1735 mutable std::vector<
1736 boost::unordered_map<SimplexId, std::list<ImplicitCluster>::iterator>>
#define TTK_TRIANGULATION_INTERNAL(NAME)
int SimplexId
Identifier type for simplices of any dimension.
std::vector< std::vector< SimplexId > > vertexStarList_
std::vector< std::array< SimplexId, 3 > > triangleList_
std::vector< std::vector< SimplexId > > triangleLinkList_
std::vector< std::vector< SimplexId > > vertexNeighborList_
std::vector< std::vector< SimplexId > > vertexLinkList_
std::vector< std::vector< SimplexId > > edgeLinkList_
std::vector< std::vector< SimplexId > > vertexEdgeList_
virtual int preconditionEdges()
std::vector< std::vector< SimplexId > > triangleEdgeVector_
std::vector< std::vector< SimplexId > > triangleStarList_
virtual SimplexId getCellTriangleNumber(const SimplexId &cellId) const
std::vector< std::vector< SimplexId > > cellNeighborList_
std::vector< std::array< SimplexId, 2 > > edgeList_
std::vector< std::vector< SimplexId > > edgeStarList_
virtual int preconditionTriangles()
std::vector< std::vector< SimplexId > > vertexTriangleList_
std::vector< std::vector< SimplexId > > cellEdgeVector_
std::vector< std::vector< SimplexId > > cellTriangleVector_
std::vector< std::vector< SimplexId > > edgeTriangleList_
int setInputCells(const SimplexId &cellNumber, const LongSimplexId *cellArray)
SimplexId findNodeIndex(SimplexId id, SIMPLEX_ID idType) const
void initCache(const float ratio=0.2)
SimplexId TTK_TRIANGULATION_INTERNAL getCellVertexNumber(const SimplexId &cellId) const override
int preconditionCellEdgesInternal() override
int preconditionVertexLinksInternal() override
int preconditionEdgesInternal() override
int preconditionBoundaryTrianglesInternal() override
int preconditionEdgeLinksInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL getNumberOfVertices() const override
const std::vector< std::vector< SimplexId > > * getCellEdgesInternal() override
int resetCache(int option)
SimplexId getVertexEdgeNumberInternal(const SimplexId &vertexId) const override
const std::vector< std::array< SimplexId, 2 > > *TTK_TRIANGULATION_INTERNAL getEdges() override
int TTK_TRIANGULATION_INTERNAL getCellNeighbor(const SimplexId &cellId, const int &localNeighborId, SimplexId &neighborId) const override
int getClusterEdgeStars(ImplicitCluster *const nodePtr) const
int getClusterVertexTriangles(ImplicitCluster *const nodePtr) const
int preconditionBoundaryVerticesInternal() override
~CompactTriangulation() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getVertexNeighbors() override
SimplexId getCellTriangleNumberInternal(const SimplexId &cellId) const override
int getTriangleVertexInternal(const SimplexId &triangleId, const int &localVertexId, SimplexId &vertexId) const override
int TTK_TRIANGULATION_INTERNAL getTriangleLink(const SimplexId &triangleId, const int &localLinkId, SimplexId &linkId) const override
SimplexId getCellEdgeNumberInternal(const SimplexId &cellId) const override
SimplexId getEdgeTriangleNumberInternal(const SimplexId &edgeId) const override
int getClusterVertexStars(ImplicitCluster *const nodePtr) const
int getClusterTriangleStars(ImplicitCluster *const nodePtr) const
SimplexId TTK_TRIANGULATION_INTERNAL getVertexLinkNumber(const SimplexId &vertexId) const override
CompactTriangulation & operator=(const CompactTriangulation &rhs)
int TTK_TRIANGULATION_INTERNAL getVertexPoint(const SimplexId &vertexId, float &x, float &y, float &z) const override
SimplexId getVertexTriangleNumberInternal(const SimplexId &vertexId) const override
const std::vector< std::vector< SimplexId > > * getCellTrianglesInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getVertexLinks() override
SimplexId TTK_TRIANGULATION_INTERNAL getTriangleStarNumber(const SimplexId &triangleId) const override
int TTK_TRIANGULATION_INTERNAL getVertexNeighbor(const SimplexId &vertexId, const int &localNeighborId, SimplexId &neighborId) const override
int preconditionCellTrianglesInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL getNumberOfCells() const override
bool TTK_TRIANGULATION_INTERNAL isVertexOnBoundary(const SimplexId &vertexId) const override
std::vector< std::vector< SimplexId > > externalCells_
SimplexId TTK_TRIANGULATION_INTERNAL getVertexNeighborNumber(const SimplexId &vertexId) const override
int TTK_TRIANGULATION_INTERNAL getEdgeStar(const SimplexId &edgeId, const int &localStarId, SimplexId &starId) const override
int getEdgeVertexInternal(const SimplexId &edgeId, const int &localVertexId, SimplexId &vertexId) const override
int getVertexTriangleInternal(const SimplexId &vertexId, const int &localTriangleId, SimplexId &triangleId) const override
std::vector< SimplexId > cellIntervals_
int countInternalTriangles(SimplexId nodeId) const
int preconditionTriangleEdgesInternal() override
SimplexId TTK_TRIANGULATION_INTERNAL getEdgeStarNumber(const SimplexId &edgeId) const override
int getClusterTriangleEdges(ImplicitCluster *const nodePtr) const
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getCellNeighbors() override
CompactTriangulation & operator=(CompactTriangulation &&rhs)=default
int getClusterCellTriangles(ImplicitCluster *const nodePtr) const
int getBoundaryCells(ImplicitCluster *const nodePtr, const SimplexId dim=2) const
SimplexId TTK_TRIANGULATION_INTERNAL getCellNeighborNumber(const SimplexId &cellId) const override
std::vector< boost::unordered_map< SimplexId, std::list< ImplicitCluster >::iterator > > cacheMaps_
int buildInternalTriangleMap(ImplicitCluster *const nodePtr, bool computeInternalTriangleList, bool computeInternalTriangleMap) const
int getClusterVertexNeighbors(ImplicitCluster *const nodePtr) const
int getClusterVertexEdges(ImplicitCluster *const nodePtr) const
const std::vector< std::vector< SimplexId > > * getTriangleEdgesInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getEdgeStars() override
std::vector< SimplexId > vertexIntervals_
int preconditionTrianglesInternal() override
int getCellTriangleInternal(const SimplexId &cellId, const int &localTriangleId, SimplexId &triangleId) const override
CompactTriangulation(CompactTriangulation &&rhs)=default
int getTriangleEdgeInternal(const SimplexId &triangleId, const int &localEdgeId, SimplexId &edgeId) const override
int preconditionEdgeStarsInternal() override
const std::vector< std::vector< SimplexId > > * getVertexTrianglesInternal() override
int getClusterEdgeLinks(ImplicitCluster *const nodePtr) const
int preconditionVertexStarsInternal() override
int buildInternalEdgeMap(ImplicitCluster *const nodePtr, bool computeInternalEdgeList, bool computeInternalEdgeMap) const
int getClusterVertexLinks(ImplicitCluster *const nodePtr) const
int TTK_TRIANGULATION_INTERNAL getDimensionality() const override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getTriangleLinks() override
int preconditionBoundaryEdgesInternal() override
int getVertexEdgeInternal(const SimplexId &vertexId, const int &localEdgeId, SimplexId &edgeId) const override
int preconditionCellNeighborsInternal() override
int TTK_TRIANGULATION_INTERNAL getVertexStar(const SimplexId &vertexId, const int &localStarId, SimplexId &starId) const override
const std::vector< std::vector< SimplexId > > * getVertexEdgesInternal() override
std::vector< SimplexId > edgeIntervals_
int getEdgeTriangleInternal(const SimplexId &edgeId, const int &localTriangleId, SimplexId &triangleId) const override
const int * vertexIndices_
int reorderVertices(std::vector< SimplexId > &vertexMap)
SimplexId TTK_TRIANGULATION_INTERNAL getEdgeLinkNumber(const SimplexId &edgeId) const override
int TTK_TRIANGULATION_INTERNAL getCellVertex(const SimplexId &cellId, const int &localVertexId, SimplexId &vertexId) const override
int preconditionTriangleLinksInternal() override
int reorderCells(const std::vector< SimplexId > &vertexMap, const SimplexId &cellNumber, const LongSimplexId *connectivity, const LongSimplexId *offset)
SimplexId countInternalEdges(SimplexId nodeId) const
int getClusterEdgeTriangles(ImplicitCluster *const nodePtr) const
int preconditionTriangleStarsInternal() override
std::shared_ptr< CellArray > cellArray_
int preconditionVertexTrianglesInternal() override
SimplexId getTriangleEdgeNumberInternal(const SimplexId &triangleId) const override
std::vector< SimplexId > triangleIntervals_
int TTK_TRIANGULATION_INTERNAL getEdgeLink(const SimplexId &edgeId, const int &localLinkId, SimplexId &linkId) const override
int buildExternalEdgeMap(ImplicitCluster *const nodePtr) const
bool TTK_TRIANGULATION_INTERNAL isEdgeOnBoundary(const SimplexId &edgeId) const override
int setInputPoints(const SimplexId &pointNumber, const void *pointSet, const int *indexArray, const bool &doublePrecision=false)
int getCellEdgeInternal(const SimplexId &cellId, const int &localEdgeId, SimplexId &edgeId) const override
int preconditionEdgeTrianglesInternal() override
int getClusterCellNeighbors(ImplicitCluster *const nodePtr) const
int preconditionVertexEdgesInternal() override
SimplexId getNumberOfTrianglesInternal() const override
int preconditionVertexNeighborsInternal() override
bool isEmpty() const override
SimplexId getNumberOfEdgesInternal() const override
SimplexId TTK_TRIANGULATION_INTERNAL getVertexStarNumber(const SimplexId &vertexId) const override
SimplexId TTK_TRIANGULATION_INTERNAL getTriangleLinkNumber(const SimplexId &triangleId) const override
int getClusterTriangleLinks(ImplicitCluster *const nodePtr) const
std::vector< std::list< ImplicitCluster > > caches_
ImplicitCluster * searchCache(const SimplexId &nodeId, const SimplexId reservedId=0) const
int buildExternalTriangleMap(ImplicitCluster *const nodePtr) const
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getEdgeLinks() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getTriangleStars() override
const std::vector< std::array< SimplexId, 3 > > *TTK_TRIANGULATION_INTERNAL getTriangles() override
bool TTK_TRIANGULATION_INTERNAL isTriangleOnBoundary(const SimplexId &triangleId) const override
int TTK_TRIANGULATION_INTERNAL getTriangleStar(const SimplexId &triangleId, const int &localStarId, SimplexId &starId) const override
const std::vector< std::vector< SimplexId > > * getEdgeTrianglesInternal() override
const std::vector< std::vector< SimplexId > > *TTK_TRIANGULATION_INTERNAL getVertexStars() override
int getClusterTetraEdges(ImplicitCluster *const nodePtr) const
int TTK_TRIANGULATION_INTERNAL getVertexLink(const SimplexId &vertexId, const int &localLinkId, SimplexId &linkId) const override
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Replacement for std::vector<std::vector<SimplexId>>
void copyTo(std::vector< std::vector< SimplexId > > &dst, int threadNumber=1) const
Copy buffers to a std::vector<std::vector<SimplexId>>
SimplexId get(SimplexId id, SimplexId local) const
Returns the data inside the sub-vectors.
SimplexId size(SimplexId id) const
Get the size of a particular sub-vector.
bool empty() const
If the underlying buffers are empty.
ImplicitCluster()=default
ImplicitCluster(SimplexId id)
friend class CompactTriangulation
~ImplicitCluster()=default
TTK base package defining the standard types.
int SimplexId
Identifier type for simplices of any dimension.
int ThreadId
Identifier type for threads (i.e. with OpenMP).
long long int LongSimplexId
Identifier type for simplices of any dimension.
T end(std::pair< T, T > &p)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)