6#include <unordered_map>
7#include <unordered_set>
9#include <hrleSparseCellIterator.hpp>
10#include <hrleSparseStarIterator.hpp>
20using namespace viennacore;
33 SmartPointer<Mesh<T>>
mesh =
nullptr;
43 return x == o.
x &&
y == o.
y &&
z == o.
z;
50 uint64_t a = (uint64_t)(uint32_t)k.
x;
51 uint64_t b = (uint64_t)(uint32_t)k.
y;
52 uint64_t c = (uint64_t)(uint32_t)k.
z;
53 uint64_t h = a * 0x9E3779B185EBCA87ULL;
54 h ^= b + 0xC2B2AE3D27D4EB4FULL + (h << 6) + (h >> 2);
55 h ^= c + 0x165667B19E3779F9ULL + (h << 6) + (h >> 2);
73 static constexpr unsigned int corner0[12] = {0, 1, 2, 0, 4, 5,
75 static constexpr unsigned int corner1[12] = {1, 3, 3, 2, 5, 7,
77 static constexpr unsigned int direction[12] = {0, 1, 0, 1, 0, 1,
85 SmartPointer<
Mesh<T>> passedMesh,
double mnd = 0.05,
105 Logger::getInstance()
106 .addError(
"No level set was passed to ToSurfaceMesh.")
110 if (
mesh ==
nullptr) {
111 Logger::getInstance()
112 .addError(
"No mesh was passed to ToSurfaceMesh.")
118 VIENNACORE_LOG_WARNING(
119 "ToSurfaceMesh: Level set is empty. No mesh will be created.");
126 VIENNACORE_LOG_WARNING(
"Levelset is less than 2 layers wide. Expanding "
127 "levelset to 2 layers.");
134 typedef std::map<hrleIndex, unsigned> nodeContainerType;
136 nodeContainerType nodes[
D];
137 nodeContainerType faceNodes[
D];
138 nodeContainerType cornerNodes;
140 typename nodeContainerType::iterator nodeIt;
144 std::vector<std::vector<unsigned>> newDataSourceIds;
147 newDataSourceIds.resize(1);
155 normalCalculator.
setMaxValue(std::numeric_limits<T>::max());
156 normalCalculator.
apply();
162 for (viennahrle::ConstSparseCellIterator<hrleDomainType> cellIt(
164 !cellIt.isFinished(); cellIt.next()) {
168 for (
int u = 0; u <
D; u++) {
169 while (!nodes[u].empty() &&
170 nodes[u].begin()->first <
hrleIndex(cellIt.getIndices()))
171 nodes[u].erase(nodes[u].begin());
173 while (!faceNodes[u].empty() &&
174 faceNodes[u].begin()->first <
hrleIndex(cellIt.getIndices()))
175 faceNodes[u].erase(faceNodes[u].begin());
178 while (!cornerNodes.empty() &&
179 cornerNodes.begin()->first <
hrleIndex(cellIt.getIndices()))
180 cornerNodes.erase(cornerNodes.begin());
187 bool hasZero =
false;
188 for (
int i = 0; i < (1 <<
D); i++) {
189 T val = cellIt.getCorner(i).getValue();
199 if (signs == (1 << (1 <<
D)) - 1 && !hasZero)
203 bool perfectCornerFound =
false;
209 for (
int i = 0; i < (1 <<
D); ++i) {
210 T val = cellIt.getCorner(i).getValue();
229 if constexpr (
D == 2) {
231 cellIt, countNeg, countPos, negMask, posMask, nodes,
232 &newDataSourceIds[0], valueIt);
233 }
else if constexpr (
D == 3) {
234 if (countNeg == 2 || countPos == 2) {
237 cellIt, countNeg, countPos, negMask, posMask, nodes, faceNodes,
238 &newDataSourceIds[0], valueIt);
239 }
else if (countNeg == 1 || countPos == 1) {
242 cellIt, countNeg, countPos, negMask, posMask, nodes,
243 cornerNodes, faceNodes, &newDataSourceIds[0], valueIt);
244 }
else if (countNeg == 3 || countPos == 3) {
247 cellIt, countNeg, countPos, negMask, posMask, nodes, faceNodes,
248 &newDataSourceIds[0], valueIt);
255 if (perfectCornerFound)
258 if constexpr (
D == 3) {
262 for (
int axis = 0; axis < 3; ++axis) {
263 for (
int d = 0; d < 2; ++d) {
268 auto it = faceNodes[axis].find(faceKey);
269 if (it != faceNodes[axis].end()) {
270 const unsigned faceNodeId = it->second;
271 const int *Triangles =
274 for (; Triangles[0] != -1; Triangles += 3) {
275 std::vector<unsigned> face_edge_nodes;
276 for (
int i = 0; i < 3; ++i) {
277 int edge = Triangles[i];
281 (((c0 >> axis) & 1) == d) && (((c1 >> axis) & 1) == d);
283 face_edge_nodes.push_back(
284 getNode(cellIt, edge, nodes, &newDataSourceIds[0]));
287 if (face_edge_nodes.size() == 2) {
289 {face_edge_nodes[0], face_edge_nodes[1], faceNodeId});
299 const int *Triangles =
303 for (; Triangles[0] != -1; Triangles +=
D) {
304 std::array<unsigned, D> nodeNumbers;
307 for (
int n = 0; n <
D; n++) {
308 const int edge = Triangles[n];
318 d += viennahrle::BitMaskToIndex<D>(p0);
320 nodeIt = nodes[dir].find(d);
321 if (nodeIt != nodes[dir].end()) {
322 nodeNumbers[n] = nodeIt->second;
324 nodeNumbers[n] =
getNode(cellIt, edge, nodes, &newDataSourceIds[0]);
336 mesh->getPointData().translateFromMultiData(
344 for (
auto &node :
mesh->nodes) {
345 for (
int i = 0; i <
D; ++i)
348 for (
int i = 0; i <
D; ++i) {
357 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
int edge) {
363 for (
int z = 0; z <
D; z++) {
365 cc[z] =
static_cast<double>(cellIt.getIndices(z) +
366 viennahrle::BitMaskToIndex<D>(p0)[z]);
368 T d0 = cellIt.getCorner(p0).getValue();
369 T d1 = cellIt.getCorner(p1).getValue();
371 cc[z] =
static_cast<T>(cellIt.getIndices(z)) +
T(0.5);
373 if (std::abs(d0) <= std::abs(d1)) {
374 cc[z] =
static_cast<T>(cellIt.getIndices(z)) + (d0 / (d0 - d1));
376 cc[z] =
static_cast<T>(cellIt.getIndices(z) + 1) - (d1 / (d1 - d0));
379 cc[z] = std::max(cc[z], cellIt.getIndices(z) +
epsilon);
380 cc[z] = std::min(cc[z], (cellIt.getIndices(z) + 1) -
epsilon);
386 unsigned getNode(viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
387 int edge, std::map<hrleIndex, unsigned> *nodes,
388 std::vector<unsigned> *newDataSourceIds) {
394 d += viennahrle::BitMaskToIndex<D>(p0);
396 if (nodes[dir].find(d) != nodes[dir].end()) {
397 return nodes[dir][d];
403 std::size_t currentPointId = 0;
404 T d0 = cellIt.getCorner(p0).getValue();
405 T d1 = cellIt.getCorner(p1).getValue();
406 if (std::abs(d0) <= std::abs(d1)) {
407 currentPointId = cellIt.getCorner(p0).getPointId();
409 currentPointId = cellIt.getCorner(p1).getPointId();
412 newDataSourceIds->push_back(currentPointId);
415 nodes[dir][d] = nodeId;
423 int axis,
bool isHighFace,
unsigned faceNodeId,
424 std::map<hrleIndex, unsigned> *nodes,
428 hrleIndex neighborIdx = cellIt.getIndices();
434 if (neighborIdx < cellIt.getIndices()) {
437 for (
int i = 0; i < 8; ++i) {
438 hrleIndex cIdx = neighborIdx + viennahrle::BitMaskToIndex<D>(i);
439 if (!grid.isOutsideOfDomain(cIdx)) {
440 valueIt.goToIndices(cIdx);
441 if (valueIt.getValue() >= 0)
446 if (nSigns != 0 && nSigns != 255) {
448 int nFaceD = isHighFace ? 0 : 1;
450 for (; nTriangles[0] != -1; nTriangles += 3) {
451 std::vector<unsigned> face_edge_nodes;
452 for (
int i = 0; i < 3; ++i) {
453 int edge = nTriangles[i];
456 bool onFace = (((c0 >> axis) & 1) == nFaceD) &&
457 (((c1 >> axis) & 1) == nFaceD);
461 hrleIndex d = neighborIdx + viennahrle::BitMaskToIndex<D>(p0);
462 auto itN = nodes[dir].find(d);
463 if (itN != nodes[dir].end()) {
464 face_edge_nodes.push_back(itN->second);
468 if (face_edge_nodes.size() == 2) {
470 face_edge_nodes[0], face_edge_nodes[1], faceNodeId});
479 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
480 int transform, Vec3D<T> &cornerPos)
const {
482 "Normal vector data must be computed for sharp corner generation.");
484 auto getTransformedGradient = [&](
int cornerID) {
485 auto corner = cellIt.getCorner(cornerID ^ transform);
486 if (corner.isDefined()) {
487 auto normal = (*normalVectorData)[corner.getPointId()];
489 if ((transform & 1) != 0)
490 normal[0] = -normal[0];
491 if ((transform & 2) != 0)
492 normal[1] = -normal[1];
498 Vec3D<T> norm1 = getTransformedGradient(1);
499 Vec3D<T> norm2 = getTransformedGradient(2);
501 if (std::abs(DotProduct(norm1, norm2)) >= 0.8) {
505 auto calculateNodePos = [&](
int edge, Vec3D<T> &pos) {
510 for (
int z = 0; z <
D; z++) {
512 pos[z] =
static_cast<T>(viennahrle::BitMaskToIndex<D>(p0)[z]);
514 const T d0 = cellIt.getCorner(p0 ^ transform).getValue();
515 const T d1 = cellIt.getCorner(p1 ^ transform).getValue();
519 if (std::abs(d0) <= std::abs(d1)) {
520 pos[z] = (d0 / (d0 - d1));
522 pos[z] =
T(1.0) - (d1 / (d1 - d0));
525 pos[z] = std::max(pos[z],
epsilon);
526 pos[z] = std::min(pos[z],
T(1.0) -
epsilon);
532 calculateNodePos(0, pX);
533 calculateNodePos(3, pY);
535 double d1 = DotProduct(norm1, pX);
536 double d2 = DotProduct(norm2, pY);
537 double det = norm1[0] * norm2[1] - norm1[1] * norm2[0];
539 if (std::abs(det) > 1e-6 && std::isfinite(det)) {
540 cornerPos[0] = (d1 * norm2[1] - d2 * norm1[1]) / det;
541 cornerPos[1] = (d2 * norm1[0] - d1 * norm2[0]) / det;
543 for (
int i = 0; i <
D; ++i)
544 cornerPos[i] = (pX[i] + pY[i]) *
T(0.5);
553 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
int countNeg,
554 int countPos,
int negMask,
int posMask,
555 std::map<hrleIndex, unsigned> *nodes,
560 for (
int i = 0; i < 4; ++i)
561 if ((negMask >> i) & 1) {
565 }
else if (countPos == 1) {
566 for (
int i = 0; i < 4; ++i)
567 if ((posMask >> i) & 1) {
573 if (cornerIdx != -1) {
575 bool isSharedCorner =
false;
577 cellIt.getIndices() + viennahrle::BitMaskToIndex<D>(cornerIdx);
578 T pVal = cellIt.getCorner(cornerIdx).getValue();
581 for (
int i = 0; i < (1 <<
D); ++i) {
586 for (
int k = 0; k <
D; ++k)
587 neighborIndices[k] = pIdx[k] - ((i >> k) & 1);
589 bool neighborIsCorner =
true;
592 for (
int k = 0; k <
D; ++k) {
593 int neighborLocal = i ^ (1 << k);
595 neighborIndices + viennahrle::BitMaskToIndex<D>(neighborLocal);
596 if (grid.isOutsideOfDomain(checkIdx)) {
597 checkIdx = grid.globalIndices2LocalIndices(checkIdx);
599 valueIt.goToIndices(checkIdx);
600 T nVal = valueIt.getValue();
601 bool pSign = (pVal >= 0);
602 bool nSign = (nVal >= 0);
603 if (pSign == nSign) {
604 neighborIsCorner =
false;
609 if (neighborIsCorner) {
610 isSharedCorner =
true;
615 if (!isSharedCorner) {
616 Vec3D<T> cornerPos{};
621 if ((cornerIdx & 1) != 0)
622 cornerPos[0] =
T(1.0) - cornerPos[0];
623 if ((cornerIdx & 2) != 0)
624 cornerPos[1] =
T(1.0) - cornerPos[1];
627 for (
int i = 0; i <
D; ++i) {
628 cornerPos[i] = (cornerPos[i] + cellIt.getIndices(i));
632 int edgeX = -1, edgeY = -1;
633 if (cornerIdx == 0) {
636 }
else if (cornerIdx == 1) {
639 }
else if (cornerIdx == 2) {
642 }
else if (cornerIdx == 3) {
647 unsigned nX =
getNode(cellIt, edgeX, nodes, newDataSourceIds);
648 unsigned nY =
getNode(cellIt, edgeY, nodes, newDataSourceIds);
655 assert(!newDataSourceIds->empty());
656 newDataSourceIds->push_back(
657 newDataSourceIds->back());
672 int countNeg,
int countPos,
int negMask,
int posMask,
673 std::map<hrleIndex, unsigned> *nodes,
674 std::map<hrleIndex, unsigned> *faceNodes,
675 std::vector<unsigned> *newDataSourceIds,
678 bool inverted =
false;
683 }
else if (countPos == 3) {
691 int C = -1, A = -1, B = -1;
692 for (
int i = 0; i < 8; ++i) {
693 if ((mask >> i) & 1) {
695 for (
int k = 0; k < 3; ++k) {
696 if ((mask >> (i ^ (1 << k))) & 1)
699 if (neighbors == 2) {
710 for (
int k = 0; k < 3; ++k) {
711 int n = C ^ (1 << k);
712 if ((mask >> n) & 1) {
721 "Normal vector data must be computed for sharp feature generation.");
723 auto getNormal = [&](
int idx) {
724 auto corner = cellIt.getCorner(idx);
725 if (corner.isDefined()) {
726 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
728 for (
int i = 0; i < 3; ++i)
737 while ((C ^ A) != (1 << axisA))
740 while ((C ^ B) != (1 << axisB))
742 int axisZ = 3 - axisA - axisB;
745 auto calculateFace = [&](
int axis,
int corner,
int n1,
746 int n2) -> std::optional<Vec3D<T>> {
748 if ((corner >> axis) & 1)
751 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
752 unsigned nodeId = faceNodes[axis][faceIdx];
753 Vec3D<T> pos =
mesh->nodes[nodeId];
754 for (
int d = 0; d < 3; ++d)
755 pos[d] = pos[d] - cellIt.getIndices(d);
759 Vec3D<T> N1 = getNormal(n1);
760 Vec3D<T> N2 = getNormal(n2);
762 if (std::abs(DotProduct(N1, N2)) >= 0.8)
765 int u = (axis + 1) % 3;
766 int v = (axis + 2) % 3;
769 while ((corner ^ n1) != (1 << axis1))
771 T t1 =
getInterp(corner, n1, cellIt, inverted);
774 while ((corner ^ n2) != (1 << axis2))
776 T t2 =
getInterp(corner, n2, cellIt, inverted);
779 P1[0] = ((corner >> 0) & 1);
780 P1[1] = ((corner >> 1) & 1);
781 P1[2] = ((corner >> 2) & 1);
782 P1[axis1] = ((corner >> axis1) & 1) ? (1.0 - t1) : t1;
785 P2[0] = ((corner >> 0) & 1);
786 P2[1] = ((corner >> 1) & 1);
787 P2[2] = ((corner >> 2) & 1);
788 P2[axis2] = ((corner >> axis2) & 1) ? (1.0 - t2) : t2;
790 double det = N1[u] * N2[v] - N1[v] * N2[u];
791 if (std::abs(det) < 1e-6)
794 double c1 = N1[u] * P1[u] + N1[v] * P1[v];
795 double c2 = N2[u] * P2[u] + N2[v] * P2[v];
799 P[u] = (c1 * N2[v] - c2 * N1[v]) / det;
800 P[v] = (c2 * N1[u] - c1 * N2[u]) / det;
802 if (Norm2(P - P1) > 9.0 || Norm2(P - P2) > 9.0)
808 auto commitFace = [&](
int axis,
int corner, Vec3D<T> P) ->
unsigned {
810 if ((corner >> axis) & 1)
813 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
814 return faceNodes[axis][faceIdx];
817 Vec3D<T> globalP = P;
818 for (
int d = 0; d < 3; ++d)
819 globalP[d] = (globalP[d] + cellIt.getIndices(d));
821 faceNodes[axis][faceIdx] = nodeId;
823 newDataSourceIds->push_back(cellIt.getCorner(corner).getPointId());
831 int D_corner = A ^ (1 << axisB);
832 int A_z = A ^ (1 << axisZ);
833 int B_z = B ^ (1 << axisZ);
834 int C_z = C ^ (1 << axisZ);
836 auto P_A_opt = calculateFace(axisA, A, D_corner, A_z);
837 auto P_B_opt = calculateFace(axisB, B, D_corner, B_z);
839 if (!P_A_opt.has_value() || !P_B_opt.has_value())
842 auto const &P_A = P_A_opt.value();
843 auto const &P_B = P_B_opt.value();
847 S[axisA] = P_B[axisA];
848 S[axisB] = P_A[axisB];
849 S[axisZ] = (P_A[axisZ] + P_B[axisZ]) * 0.5;
853 for (
int i = 0; i < 3; ++i)
854 gS[i] = (gS[i] + cellIt.getIndices(i));
858 newDataSourceIds->push_back(cellIt.getCorner(C).getPointId());
862 S_Face[axisZ] =
static_cast<T>((C >> axisZ) & 1);
865 auto faceIdxZ = cellIt.getIndices();
866 if ((C >> axisZ) & 1)
869 if (faceNodes[axisZ].find(faceIdxZ) != faceNodes[axisZ].end()) {
870 nS_Face = faceNodes[axisZ][faceIdxZ];
873 Vec3D<T> gS = S_Face;
874 for (
int i = 0; i < 3; ++i)
875 gS[i] = (gS[i] + cellIt.getIndices(i));
878 faceNodes[axisZ][faceIdxZ] = nS_Face;
880 newDataSourceIds[0].push_back(cellIt.getCorner(C).getPointId());
887 unsigned nNA = commitFace(axisA, A, P_A);
888 unsigned nNB = commitFace(axisB, B, P_B);
891 auto getEdgeNode = [&](
int v1,
int v2) {
893 for (
int e = 0; e < 12; ++e) {
894 if ((
corner0[e] ==
static_cast<unsigned>(v1) &&
895 corner1[e] ==
static_cast<unsigned>(v2)) ||
896 (
corner0[e] ==
static_cast<unsigned>(v2) &&
897 corner1[e] ==
static_cast<unsigned>(v1))) {
902 return this->
getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
905 unsigned nI_AD = getEdgeNode(A, D_corner);
906 unsigned nI_AZ = getEdgeNode(A, A_z);
907 unsigned nI_BD = getEdgeNode(B, D_corner);
908 unsigned nI_BZ = getEdgeNode(B, B_z);
909 unsigned nI_CZ = getEdgeNode(C, C_z);
912 int parity = (C & 1) + ((C >> 1) & 1) + ((C >> 2) & 1);
913 bool flip = (parity % 2 != 0) ^ inverted;
916 while ((A ^ C) != (1 << vA))
919 while ((B ^ C) != (1 << vB))
921 bool is_cyclic = ((vA + 1) % 3 == vB);
924 std::swap(nI_AD, nI_BD);
925 std::swap(nI_AZ, nI_BZ);
954 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
955 int transform,
int axis,
bool inverted,
956 std::map<hrleIndex, unsigned> *nodes,
957 std::map<hrleIndex, unsigned> *faceNodes,
960 "Normal vector data must be computed for sharp edge generation.");
963 auto toCanonical = [&](Vec3D<T> v) {
979 res = Vec3D<T>{v[1], v[2], v[0]};
981 res = Vec3D<T>{v[2], v[0], v[1]};
986 auto fromCanonical = [&](Vec3D<T> v) {
992 res = Vec3D<T>{v[2], v[0], v[1]};
994 res = Vec3D<T>{v[1], v[2], v[0]};
997 res[0] =
T(1.0) - res[0];
999 res[1] =
T(1.0) - res[1];
1001 res[2] =
T(1.0) - res[2];
1005 auto getNormal = [&](
int idx) {
1006 auto corner = cellIt.getCorner(idx);
1007 if (corner.isDefined()) {
1008 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1010 for (
int i = 0; i < 3; ++i)
1020 unsigned faceNodeIds[2];
1021 bool nodeExists[2] = {
false,
false};
1024 auto mapIdx = [&](
int c) {
1027 int x = (c & 1), y = (c >> 1) & 1, z = (c >> 2) & 1;
1033 }
else if (axis == 1) {
1051 return gx | (gy << 1) | (gz << 2);
1055 for (
int k = 0; k < 2; ++k) {
1062 hrleIndex faceIdx = cellIt.getIndices();
1066 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1070 auto it = faceNodes[axis].find(faceIdx);
1071 if (it != faceNodes[axis].end()) {
1072 nodeExists[k] =
true;
1073 faceNodeIds[k] = it->second;
1074 Vec3D<T> relPos =
mesh->nodes[faceNodeIds[k]];
1075 for (
int d = 0; d < 3; ++d) {
1076 relPos[d] -=
static_cast<T>(cellIt.getIndices(d));
1078 P[k] = toCanonical(relPos);
1094 Vec3D<T> n_y = toCanonical(getNormal(c_y));
1095 Vec3D<T> n_z = toCanonical(getNormal(c_z));
1097 if (std::abs(DotProduct(n_y, n_z)) >= 0.8) {
1105 T t_y =
getInterp(c0, c_y, cellIt, inverted);
1106 T t_z =
getInterp(c0, c_z, cellIt, inverted);
1116 double det = n_y[1] * n_z[2] - n_y[2] * n_z[1];
1117 if (std::abs(det) < 1e-6)
1120 double d1 = n_y[1] * t_y;
1121 double d2 = n_z[2] * t_z;
1123 double y = (d1 * n_z[2] - d2 * n_y[2]) / det;
1124 double z = (d2 * n_y[1] - d1 * n_z[1]) / det;
1126 P[k] = Vec3D<T>{(
T)k, (
T)y, (
T)z};
1130 if (Norm2(P[k] - Vec3D<T>{(
T)k, t_y, 0}) > 9.0 ||
1131 Norm2(P[k] - Vec3D<T>{(
T)k, 0, t_z}) > 9.0) {
1138 for (
int k = 0; k < 2; ++k) {
1139 if (!nodeExists[k]) {
1141 hrleIndex faceIdx = cellIt.getIndices();
1142 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1146 Vec3D<T> globalP = fromCanonical(P[k]);
1147 for (
int d = 0; d < 3; ++d)
1148 globalP[d] = (globalP[d] + cellIt.getIndices(d));
1151 faceNodes[axis][faceIdx] = faceNodeIds[k];
1153 newDataSourceIds->push_back(
1154 cellIt.getCorner(mapIdx(c0)).getPointId());
1169 auto getEdgeNode = [&](
int c1,
int c2) {
1171 int g1 = mapIdx(c1);
1172 int g2 = mapIdx(c2);
1176 for (
int e = 0; e < 12; ++e) {
1177 if ((
corner0[e] ==
static_cast<unsigned>(g1) &&
1178 corner1[e] ==
static_cast<unsigned>(g2)) ||
1179 (
corner0[e] ==
static_cast<unsigned>(g2) &&
1180 corner1[e] ==
static_cast<unsigned>(g1))) {
1185 return this->
getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1188 unsigned n02 = getEdgeNode(0, 2);
1189 unsigned n04 = getEdgeNode(0, 4);
1190 unsigned n13 = getEdgeNode(1, 3);
1191 unsigned n15 = getEdgeNode(1, 5);
1210 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
int countNeg,
1211 int countPos,
int negMask,
int posMask,
1212 std::map<hrleIndex, unsigned> *nodes,
1213 std::map<hrleIndex, unsigned> *faceNodes,
1216 bool inverted =
false;
1218 if (countNeg == 2) {
1221 }
else if (countPos == 2) {
1229 int v1 = -1, v2 = -1;
1230 for (
int i = 0; i < 8; ++i) {
1231 if ((mask >> i) & 1) {
1240 if ((diff & (diff - 1)) != 0)
1244 while ((diff >> (axis + 1)) > 0)
1251 nodes, faceNodes, newDataSourceIds,
1258 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
1259 int transform,
bool inverted, std::map<hrleIndex, unsigned> *nodes,
1260 std::map<hrleIndex, unsigned> *faceNodes,
1262 Vec3D<T> &cornerPos) {
1268 "Normal vector data must be available for sharp corner generation");
1271 auto toCanonical = [&](Vec3D<T> v) {
1282 auto fromCanonical = [&](Vec3D<T> v) {
1284 v[0] =
T(1.0) - v[0];
1286 v[1] =
T(1.0) - v[1];
1288 v[2] =
T(1.0) - v[2];
1292 auto getNormal = [&](
int idx) {
1293 auto corner = cellIt.getCorner(idx);
1294 if (corner.isDefined()) {
1295 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1297 for (
int i = 0; i < 3; ++i)
1305 Vec3D<T> norm1, norm2, norm3;
1307 Vec3D<T> P1, P2, P4;
1309 const int g1 = transform ^ 1;
1310 const int g2 = transform ^ 2;
1311 const int g4 = transform ^ 4;
1313 norm1 = toCanonical(getNormal(g1));
1314 norm2 = toCanonical(getNormal(g2));
1315 norm3 = toCanonical(getNormal(g4));
1317 const T t1 =
getInterp(transform, g1, cellIt, inverted);
1318 const T t2 =
getInterp(transform, g2, cellIt, inverted);
1319 const T t4 =
getInterp(transform, g4, cellIt, inverted);
1325 d1 = DotProduct(norm1, P1);
1326 d2 = DotProduct(norm2, P2);
1327 d3 = DotProduct(norm3, P4);
1329 if (std::abs(DotProduct(norm1, norm2)) >= 0.8 ||
1330 std::abs(DotProduct(norm1, norm3)) >= 0.8 ||
1331 std::abs(DotProduct(norm2, norm3)) >= 0.8) {
1344 const Vec3D<T> c23 = CrossProduct(norm2, norm3);
1345 const Vec3D<T> c31 = CrossProduct(norm3, norm1);
1346 const Vec3D<T> c12 = CrossProduct(norm1, norm2);
1348 const double det = DotProduct(norm1, c23);
1353 const T invDet = 1.0 / det;
1355 for (
int i = 0; i < 3; ++i) {
1356 S[i] = (d1 * c23[i] + d2 * c31[i] + d3 * c12[i]) * invDet;
1359 if (Norm2(S - P1) > 9.0 || Norm2(S - P2) > 9.0 || Norm2(S - P4) > 9.0) {
1364 Vec3D<T> globalS = fromCanonical(S);
1365 for (
int d = 0; d < 3; ++d) {
1366 cornerPos[d] = (globalS[d] + cellIt.getIndices(d));
1375 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
int countNeg,
1376 int countPos,
int negMask,
int posMask,
1377 std::map<hrleIndex, unsigned> *nodes,
1378 std::map<hrleIndex, unsigned> &cornerNodes,
1379 std::map<hrleIndex, unsigned> *faceNodes,
1382 bool inverted =
false;
1386 if (countNeg == 1) {
1388 while (!((negMask >> transform) & 1))
1391 }
else if (countPos == 1) {
1393 while (!((posMask >> transform) & 1))
1398 if (transform == -1)
1403 faceNodes, newDataSourceIds, valueIt,
1406 hrleIndex cornerIdx = cellIt.getIndices();
1407 cornerIdx += viennahrle::BitMaskToIndex<D>(transform);
1409 auto it = cornerNodes.find(cornerIdx);
1410 if (it != cornerNodes.end()) {
1411 nCorner = it->second;
1414 cornerNodes[cornerIdx] = nCorner;
1420 newDataSourceIds->push_back(cellIt.getCorner(transform).getPointId());
1424 auto getEdgeNode = [&](
int neighbor) {
1428 for (
int e = 0; e < 12; ++e) {
1429 if ((
corner0[e] ==
static_cast<unsigned>(v1) &&
1430 corner1[e] ==
static_cast<unsigned>(v2)) ||
1431 (
corner0[e] ==
static_cast<unsigned>(v2) &&
1432 corner1[e] ==
static_cast<unsigned>(v1))) {
1437 return this->
getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1440 unsigned n1 = getEdgeNode(transform ^ 1);
1441 unsigned n2 = getEdgeNode(transform ^ 2);
1442 unsigned n4 = getEdgeNode(transform ^ 4);
1447 auto toCanonical = [&](Vec3D<T> v) {
1457 auto fromCanonical = [&](Vec3D<T> v) {
1459 v[0] =
T(1.0) - v[0];
1461 v[1] =
T(1.0) - v[1];
1463 v[2] =
T(1.0) - v[2];
1467 auto getNormal = [&](
int idx) {
1468 auto corner = cellIt.getCorner(idx);
1469 if (corner.isDefined()) {
1470 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1472 for (
int i = 0; i < 3; ++i)
1479 Vec3D<T> norm1, norm2, norm3;
1483 int n1_idx = -1, n2_idx = -1, n3_idx = -1;
1484 int neighbors[] = {transform ^ 1, transform ^ 2, transform ^ 4};
1486 for (
int n_idx : neighbors) {
1487 T val = cellIt.getCorner(n_idx).getValue();
1489 if (val * cellIt.getCorner(transform).getValue() < 0) {
1500 norm1 = toCanonical(getNormal(transform ^ 1));
1501 norm2 = toCanonical(getNormal(transform ^ 2));
1502 norm3 = toCanonical(getNormal(transform ^ 4));
1504 d1 = norm1[0] *
getInterp(transform, transform ^ 1, cellIt, inverted);
1505 d2 = norm2[1] *
getInterp(transform, transform ^ 2, cellIt, inverted);
1506 d3 = norm3[2] *
getInterp(transform, transform ^ 4, cellIt, inverted);
1508 auto solve2x2 = [&](
double a1,
double b1,
double c1,
double a2,
1509 double b2,
double c2,
T &res1,
T &res2) {
1510 double det = a1 * b2 - a2 * b1;
1511 if (std::abs(det) < 1e-6)
1513 res1 = (c1 * b2 - c2 * b1) / det;
1514 res2 = (a1 * c2 - a2 * c1) / det;
1518 Vec3D<T> Fx, Fy, Fz;
1524 valid &= solve2x2(norm2[1], norm2[2], d2 - norm2[0] * Fx[0], norm3[1],
1525 norm3[2], d3 - norm3[0] * Fx[0], Fx[1], Fx[2]);
1526 valid &= solve2x2(norm1[0], norm1[2], d1 - norm1[1] * Fy[1], norm3[0],
1527 norm3[2], d3 - norm3[1] * Fy[1], Fy[0], Fy[2]);
1528 valid &= solve2x2(norm1[0], norm1[1], d1 - norm1[2] * Fz[2], norm2[0],
1529 norm2[1], d2 - norm2[2] * Fz[2], Fz[0], Fz[1]);
1531 auto checkBounds = [&](Vec3D<T> p) {
1532 return p[0] >= -1e-4 && p[0] <= 1.0 + 1e-4 && p[1] >= -1e-4 &&
1533 p[1] <= 1.0 + 1e-4 && p[2] >= -1e-4 && p[2] <= 1.0 + 1e-4;
1536 if (valid && checkBounds(Fx) && checkBounds(Fy) && checkBounds(Fz)) {
1537 auto getFaceNode = [&](
int axis, Vec3D<T> pos) {
1538 hrleIndex faceIdx = cellIt.getIndices();
1539 if ((transform >> axis) & 1)
1542 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
1543 return faceNodes[axis][faceIdx];
1545 Vec3D<T> globalPos = fromCanonical(pos);
1546 for (
int d = 0; d < 3; ++d)
1547 globalPos[d] = (globalPos[d] + cellIt.getIndices(d));
1549 faceNodes[axis][faceIdx] = nodeId;
1551 newDataSourceIds->push_back(
1552 cellIt.getCorner(transform).getPointId());
1558 unsigned nFx = getFaceNode(0, Fx);
1559 unsigned nFy = getFaceNode(1, Fy);
1560 unsigned nFz = getFaceNode(2, Fz);
1564 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1567 bool flip = (parity % 2 != 0) ^ inverted;
1595 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1598 bool flip = (parity % 2 != 0) ^ inverted;
1618 if constexpr (
D == 3) {
1619 return nodeNumbers[0] == nodeNumbers[1] ||
1620 nodeNumbers[0] == nodeNumbers[2] ||
1621 nodeNumbers[1] == nodeNumbers[2];
1623 return nodeNumbers[0] == nodeNumbers[1];
1628 const Vec3D<T> &nodeB,
1629 const Vec3D<T> &nodeC)
noexcept {
1630 return CrossProduct(nodeB - nodeA, nodeC - nodeA);
1634 const viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
1635 bool inverted)
const {
1636 auto getValue = [&](
int idx) {
1637 T val = cellIt.getCorner(idx).getValue();
1638 return inverted ? -val : val;
1640 const T v_a = getValue(p_a);
1641 const T v_b = getValue(p_b);
1642 if (std::abs(v_a - v_b) <
epsilon)
1644 return v_a / (v_a - v_b);
1651 auto elementI3 = [&](
const std::array<unsigned, D> &element) ->
I3 {
1652 if constexpr (
D == 2)
1653 return {(int)element[0], (
int)element[1], 0};
1655 return {(int)element[0], (
int)element[1], (int)element[2]};
1660 if constexpr (
D == 2) {
1662 -(
mesh->nodes[nodeNumbers[1]][1] -
mesh->nodes[nodeNumbers[0]][1]),
1663 mesh->nodes[nodeNumbers[1]][0] -
mesh->nodes[nodeNumbers[0]][0],
1667 mesh->nodes[nodeNumbers[1]],
1668 mesh->nodes[nodeNumbers[2]]);
1672 normal[0] * normal[0] + normal[1] * normal[1] + normal[2] * normal[2];
1674 mesh->insertNextElement(nodeNumbers);
1675 T invn =
static_cast<T>(1.) / std::sqrt(
static_cast<T>(n2));
1676 for (
int d = 0; d <
D; d++) {
1686 auto quantize = [&](
const Vec3D<T> &p) ->
I3 {
1688 return {(int)std::llround(p[0] * inv), (int)std::llround(p[1] * inv),
1689 (int)std::llround(p[2] * inv)};
1694 auto q = quantize(pos);
1697 nodeIdx = it->second;
1701 for (
int i = 0; i < 3; ++i) {
1702 mesh->nodes[nodeIdx][i] = (
mesh->nodes[nodeIdx][i] + pos[i]) *
T(0.5);
1706 unsigned newNodeId =
mesh->insertNextNode(pos);
1710 for (
int a = 0; a <
D; a++) {
1711 if (pos[a] <
mesh->minimumExtent[a])
1712 mesh->minimumExtent[a] = pos[a];
1713 if (pos[a] >
mesh->maximumExtent[a])
1714 mesh->maximumExtent[a] = pos[a];
constexpr int D
Definition Epitaxy.cpp:11
double T
Definition Epitaxy.cpp:12
static const int * polygonize2d(unsigned int signs)
signs = signs of the corners in lexicographic order (1bit per corner)
Definition lsMarchingCubes.hpp:312
static const int * polygonize3d(unsigned int signs)
signs = signs of the corners in lexicographic order (1bit per corner)
Definition lsMarchingCubes.hpp:324
This algorithm is used to compute the normal vectors for all points with level set values <= maxValue...
Definition lsCalculateNormalVectors.hpp:40
static constexpr char normalVectorsLabel[]
Definition lsCalculateNormalVectors.hpp:52
void apply()
Definition lsCalculateNormalVectors.hpp:95
void setMethod(NormalCalculationMethodEnum passedMethod)
Definition lsCalculateNormalVectors.hpp:79
void setMaxValue(const T passedMaxValue)
Definition lsCalculateNormalVectors.hpp:67
Class containing all information about the level set, including the dimensions of the domain,...
Definition lsDomain.hpp:28
viennahrle::Domain< T, D > DomainType
Definition lsDomain.hpp:33
Expands the levelSet to the specified number of layers. The largest value in the levelset is thus wid...
Definition lsExpand.hpp:17
void apply()
Apply the expansion to the specified width.
Definition lsExpand.hpp:44
This class holds an explicit mesh, which is always given in 3 dimensions. If it describes a 2D mesh,...
Definition lsMesh.hpp:22
std::vector< Vec3D< T > > VectorDataType
Definition lsPointData.hpp:25
SmartPointer< Mesh< T > > mesh
Definition lsToSurfaceMesh.hpp:33
ToSurfaceMesh(const SmartPointer< lsDomainType > passedLevelSet, SmartPointer< Mesh< T > > passedMesh, double mnd=0.05, double eps=1e-12)
Definition lsToSurfaceMesh.hpp:84
static constexpr unsigned int direction[12]
Definition lsToSurfaceMesh.hpp:77
bool generateSharpEdge3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:1209
void setMesh(SmartPointer< Mesh< T > > passedMesh)
Definition lsToSurfaceMesh.hpp:95
SmartPointer< lsDomainType > currentLevelSet
Definition lsToSurfaceMesh.hpp:67
bool generateSharpCorner2D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:552
bool generateSharpCorners
Definition lsToSurfaceMesh.hpp:38
bool generateCanonicalSharpEdge3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int transform, int axis, bool inverted, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:953
std::vector< SmartPointer< lsDomainType > > levelSets
Definition lsToSurfaceMesh.hpp:32
typename lsDomainType::DomainType hrleDomainType
Definition lsToSurfaceMesh.hpp:28
static bool triangleMisformed(const std::array< unsigned, D > &nodeNumbers) noexcept
Definition lsToSurfaceMesh.hpp:1617
std::unordered_map< I3, unsigned, I3Hash > nodeIdByBin
Definition lsToSurfaceMesh.hpp:61
const T epsilon
Definition lsToSurfaceMesh.hpp:35
viennals::Domain< T, D > lsDomainType
Definition lsToSurfaceMesh.hpp:27
viennahrle::ConstSparseIterator< hrleDomainType > ConstSparseIterator
Definition lsToSurfaceMesh.hpp:30
Vec3D< T > computeNodePosition(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int edge)
Definition lsToSurfaceMesh.hpp:356
bool generateCanonicalSharpCorner2D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int transform, Vec3D< T > &cornerPos) const
Definition lsToSurfaceMesh.hpp:478
bool generateSharpCorner3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > &cornerNodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:1374
void stitchToNeighbor(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int axis, bool isHighFace, unsigned faceNodeId, std::map< hrleIndex, unsigned > *nodes, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:422
std::unordered_set< I3, I3Hash > uniqueElements
Definition lsToSurfaceMesh.hpp:62
std::vector< std::pair< unsigned, Vec3D< T > > > matSharpCornerNodes
Definition lsToSurfaceMesh.hpp:71
static constexpr unsigned int corner1[12]
Definition lsToSurfaceMesh.hpp:75
static Vec3D< T > calculateNormal(const Vec3D< T > &nodeA, const Vec3D< T > &nodeB, const Vec3D< T > &nodeC) noexcept
Definition lsToSurfaceMesh.hpp:1627
const double minNodeDistanceFactor
Definition lsToSurfaceMesh.hpp:36
void insertElement(const std::array< unsigned, D > &nodeNumbers)
Definition lsToSurfaceMesh.hpp:1647
void setUpdatePointData(bool update)
Definition lsToSurfaceMesh.hpp:97
ToSurfaceMesh(double mnd=0.05, double eps=1e-12)
Definition lsToSurfaceMesh.hpp:81
double currentGridDelta
Definition lsToSurfaceMesh.hpp:65
unsigned insertNode(Vec3D< T > const &pos)
Definition lsToSurfaceMesh.hpp:1685
viennahrle::Index< D > hrleIndex
Definition lsToSurfaceMesh.hpp:29
bool updatePointData
Definition lsToSurfaceMesh.hpp:37
bool generateCanonicalSharpCorner3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int transform, bool inverted, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt, Vec3D< T > &cornerPos)
Definition lsToSurfaceMesh.hpp:1257
T getInterp(int p_a, int p_b, const viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, bool inverted) const
Definition lsToSurfaceMesh.hpp:1633
bool generateSharpL3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:671
PointData< T >::VectorDataType * normalVectorData
Definition lsToSurfaceMesh.hpp:68
std::vector< Vec3D< T > > currentNormals
Definition lsToSurfaceMesh.hpp:63
void scaleMesh()
Definition lsToSurfaceMesh.hpp:342
std::vector< T > currentMaterials
Definition lsToSurfaceMesh.hpp:64
void setSharpCorners(bool check)
Definition lsToSurfaceMesh.hpp:99
virtual void apply()
Definition lsToSurfaceMesh.hpp:101
T currentMaterialId
Definition lsToSurfaceMesh.hpp:66
void setLevelSet(SmartPointer< lsDomainType > passedlsDomain)
Definition lsToSurfaceMesh.hpp:90
unsigned getNode(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int edge, std::map< hrleIndex, unsigned > *nodes, std::vector< unsigned > *newDataSourceIds)
Definition lsToSurfaceMesh.hpp:386
static constexpr unsigned int corner0[12]
Definition lsToSurfaceMesh.hpp:73
#define PRECOMPILE_PRECISION_DIMENSION(className)
Definition lsPreCompileMacros.hpp:24
Definition lsAdvect.hpp:41
@ ONE_SIDED_MIN_MOD
Definition lsCalculateNormalVectors.hpp:22
Definition lsToSurfaceMesh.hpp:47
size_t operator()(const I3 &k) const
Definition lsToSurfaceMesh.hpp:48
Definition lsToSurfaceMesh.hpp:40
bool operator==(const I3 &o) const
Definition lsToSurfaceMesh.hpp:42
int x
Definition lsToSurfaceMesh.hpp:41
int y
Definition lsToSurfaceMesh.hpp:41
int z
Definition lsToSurfaceMesh.hpp:41