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;
421 template <
int Dim = D>
422 std::enable_if_t<Dim == 3, void>
424 int axis,
bool isHighFace,
unsigned faceNodeId,
425 std::map<hrleIndex, unsigned> *nodes,
429 hrleIndex neighborIdx = cellIt.getIndices();
435 if (neighborIdx < cellIt.getIndices()) {
438 for (
int i = 0; i < 8; ++i) {
439 hrleIndex cIdx = neighborIdx + viennahrle::BitMaskToIndex<D>(i);
440 if (!grid.isOutsideOfDomain(cIdx)) {
441 valueIt.goToIndices(cIdx);
442 if (valueIt.getValue() >= 0)
447 if (nSigns != 0 && nSigns != 255) {
449 int nFaceD = isHighFace ? 0 : 1;
451 for (; nTriangles[0] != -1; nTriangles += 3) {
452 std::vector<unsigned> face_edge_nodes;
453 for (
int i = 0; i < 3; ++i) {
454 int edge = nTriangles[i];
457 bool onFace = (((c0 >> axis) & 1) == nFaceD) &&
458 (((c1 >> axis) & 1) == nFaceD);
462 hrleIndex d = neighborIdx + viennahrle::BitMaskToIndex<D>(p0);
463 auto itN = nodes[dir].find(d);
464 if (itN != nodes[dir].end()) {
465 face_edge_nodes.push_back(itN->second);
469 if (face_edge_nodes.size() == 2) {
471 face_edge_nodes[0], face_edge_nodes[1], faceNodeId});
479 template <
int Dim = D>
481 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
482 int transform, Vec3D<T> &cornerPos)
const {
484 "Normal vector data must be computed for sharp corner generation.");
486 auto getTransformedGradient = [&](
int cornerID) {
487 auto corner = cellIt.getCorner(cornerID ^ transform);
488 if (corner.isDefined()) {
489 auto normal = (*normalVectorData)[corner.getPointId()];
491 if ((transform & 1) != 0)
492 normal[0] = -normal[0];
493 if ((transform & 2) != 0)
494 normal[1] = -normal[1];
500 Vec3D<T> norm1 = getTransformedGradient(1);
501 Vec3D<T> norm2 = getTransformedGradient(2);
503 if (std::abs(DotProduct(norm1, norm2)) >= 0.8) {
507 auto calculateNodePos = [&](
int edge, Vec3D<T> &pos) {
512 for (
int z = 0; z <
D; z++) {
514 pos[z] =
static_cast<T>(viennahrle::BitMaskToIndex<D>(p0)[z]);
516 const T d0 = cellIt.getCorner(p0 ^ transform).getValue();
517 const T d1 = cellIt.getCorner(p1 ^ transform).getValue();
521 if (std::abs(d0) <= std::abs(d1)) {
522 pos[z] = (d0 / (d0 - d1));
524 pos[z] =
T(1.0) - (d1 / (d1 - d0));
527 pos[z] = std::max(pos[z],
epsilon);
528 pos[z] = std::min(pos[z],
T(1.0) -
epsilon);
534 calculateNodePos(0, pX);
535 calculateNodePos(3, pY);
537 double d1 = DotProduct(norm1, pX);
538 double d2 = DotProduct(norm2, pY);
539 double det = norm1[0] * norm2[1] - norm1[1] * norm2[0];
541 if (std::abs(det) > 1e-6 && std::isfinite(det)) {
542 cornerPos[0] = (d1 * norm2[1] - d2 * norm1[1]) / det;
543 cornerPos[1] = (d2 * norm1[0] - d1 * norm2[0]) / det;
545 for (
int i = 0; i <
D; ++i)
546 cornerPos[i] = (pX[i] + pY[i]) *
T(0.5);
554 template <
int Dim = D>
556 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
int countNeg,
557 int countPos,
int negMask,
int posMask,
558 std::map<hrleIndex, unsigned> *nodes,
563 for (
int i = 0; i < 4; ++i)
564 if ((negMask >> i) & 1) {
568 }
else if (countPos == 1) {
569 for (
int i = 0; i < 4; ++i)
570 if ((posMask >> i) & 1) {
576 if (cornerIdx != -1) {
578 bool isSharedCorner =
false;
580 cellIt.getIndices() + viennahrle::BitMaskToIndex<D>(cornerIdx);
581 T pVal = cellIt.getCorner(cornerIdx).getValue();
584 for (
int i = 0; i < (1 <<
D); ++i) {
589 for (
int k = 0; k <
D; ++k)
590 neighborIndices[k] = pIdx[k] - ((i >> k) & 1);
592 bool neighborIsCorner =
true;
595 for (
int k = 0; k <
D; ++k) {
596 int neighborLocal = i ^ (1 << k);
598 neighborIndices + viennahrle::BitMaskToIndex<D>(neighborLocal);
599 if (grid.isOutsideOfDomain(checkIdx)) {
600 checkIdx = grid.globalIndices2LocalIndices(checkIdx);
602 valueIt.goToIndices(checkIdx);
603 T nVal = valueIt.getValue();
604 bool pSign = (pVal >= 0);
605 bool nSign = (nVal >= 0);
606 if (pSign == nSign) {
607 neighborIsCorner =
false;
612 if (neighborIsCorner) {
613 isSharedCorner =
true;
618 if (!isSharedCorner) {
619 Vec3D<T> cornerPos{};
624 if ((cornerIdx & 1) != 0)
625 cornerPos[0] =
T(1.0) - cornerPos[0];
626 if ((cornerIdx & 2) != 0)
627 cornerPos[1] =
T(1.0) - cornerPos[1];
630 for (
int i = 0; i <
D; ++i) {
631 cornerPos[i] = (cornerPos[i] + cellIt.getIndices(i));
635 int edgeX = -1, edgeY = -1;
636 if (cornerIdx == 0) {
639 }
else if (cornerIdx == 1) {
642 }
else if (cornerIdx == 2) {
645 }
else if (cornerIdx == 3) {
650 unsigned nX =
getNode(cellIt, edgeX, nodes, newDataSourceIds);
651 unsigned nY =
getNode(cellIt, edgeY, nodes, newDataSourceIds);
658 assert(!newDataSourceIds->empty());
659 newDataSourceIds->push_back(
660 newDataSourceIds->back());
673 template <
int Dim = D>
674 std::enable_if_t<Dim == 3, bool>
676 int countNeg,
int countPos,
int negMask,
int posMask,
677 std::map<hrleIndex, unsigned> *nodes,
678 std::map<hrleIndex, unsigned> *faceNodes,
679 std::vector<unsigned> *newDataSourceIds,
682 bool inverted =
false;
687 }
else if (countPos == 3) {
695 int C = -1, A = -1, B = -1;
696 for (
int i = 0; i < 8; ++i) {
697 if ((mask >> i) & 1) {
699 for (
int k = 0; k < 3; ++k) {
700 if ((mask >> (i ^ (1 << k))) & 1)
703 if (neighbors == 2) {
714 for (
int k = 0; k < 3; ++k) {
715 int n = C ^ (1 << k);
716 if ((mask >> n) & 1) {
725 "Normal vector data must be computed for sharp feature generation.");
727 auto getNormal = [&](
int idx) {
728 auto corner = cellIt.getCorner(idx);
729 if (corner.isDefined()) {
730 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
732 for (
int i = 0; i < 3; ++i)
741 while ((C ^ A) != (1 << axisA))
744 while ((C ^ B) != (1 << axisB))
746 int axisZ = 3 - axisA - axisB;
749 auto calculateFace = [&](
int axis,
int corner,
int n1,
750 int n2) -> std::optional<Vec3D<T>> {
752 if ((corner >> axis) & 1)
755 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
756 unsigned nodeId = faceNodes[axis][faceIdx];
757 Vec3D<T> pos =
mesh->nodes[nodeId];
758 for (
int d = 0; d < 3; ++d)
759 pos[d] = pos[d] - cellIt.getIndices(d);
763 Vec3D<T> N1 = getNormal(n1);
764 Vec3D<T> N2 = getNormal(n2);
766 if (std::abs(DotProduct(N1, N2)) >= 0.8)
769 int u = (axis + 1) % 3;
770 int v = (axis + 2) % 3;
773 while ((corner ^ n1) != (1 << axis1))
775 T t1 =
getInterp(corner, n1, cellIt, inverted);
778 while ((corner ^ n2) != (1 << axis2))
780 T t2 =
getInterp(corner, n2, cellIt, inverted);
783 P1[0] = ((corner >> 0) & 1);
784 P1[1] = ((corner >> 1) & 1);
785 P1[2] = ((corner >> 2) & 1);
786 P1[axis1] = ((corner >> axis1) & 1) ? (1.0 - t1) : t1;
789 P2[0] = ((corner >> 0) & 1);
790 P2[1] = ((corner >> 1) & 1);
791 P2[2] = ((corner >> 2) & 1);
792 P2[axis2] = ((corner >> axis2) & 1) ? (1.0 - t2) : t2;
794 double det = N1[u] * N2[v] - N1[v] * N2[u];
795 if (std::abs(det) < 1e-6)
798 double c1 = N1[u] * P1[u] + N1[v] * P1[v];
799 double c2 = N2[u] * P2[u] + N2[v] * P2[v];
803 P[u] = (c1 * N2[v] - c2 * N1[v]) / det;
804 P[v] = (c2 * N1[u] - c1 * N2[u]) / det;
806 if (Norm2(P - P1) > 9.0 || Norm2(P - P2) > 9.0)
812 auto commitFace = [&](
int axis,
int corner, Vec3D<T> P) ->
unsigned {
814 if ((corner >> axis) & 1)
817 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
818 return faceNodes[axis][faceIdx];
821 Vec3D<T> globalP = P;
822 for (
int d = 0; d < 3; ++d)
823 globalP[d] = (globalP[d] + cellIt.getIndices(d));
825 faceNodes[axis][faceIdx] = nodeId;
827 newDataSourceIds->push_back(cellIt.getCorner(corner).getPointId());
835 int D_corner = A ^ (1 << axisB);
836 int A_z = A ^ (1 << axisZ);
837 int B_z = B ^ (1 << axisZ);
838 int C_z = C ^ (1 << axisZ);
840 auto P_A_opt = calculateFace(axisA, A, D_corner, A_z);
841 auto P_B_opt = calculateFace(axisB, B, D_corner, B_z);
843 if (!P_A_opt.has_value() || !P_B_opt.has_value())
846 auto const &P_A = P_A_opt.value();
847 auto const &P_B = P_B_opt.value();
851 S[axisA] = P_B[axisA];
852 S[axisB] = P_A[axisB];
853 S[axisZ] = (P_A[axisZ] + P_B[axisZ]) * 0.5;
857 for (
int i = 0; i < 3; ++i)
858 gS[i] = (gS[i] + cellIt.getIndices(i));
862 newDataSourceIds->push_back(cellIt.getCorner(C).getPointId());
866 S_Face[axisZ] =
static_cast<T>((C >> axisZ) & 1);
869 auto faceIdxZ = cellIt.getIndices();
870 if ((C >> axisZ) & 1)
873 if (faceNodes[axisZ].find(faceIdxZ) != faceNodes[axisZ].end()) {
874 nS_Face = faceNodes[axisZ][faceIdxZ];
877 Vec3D<T> gS = S_Face;
878 for (
int i = 0; i < 3; ++i)
879 gS[i] = (gS[i] + cellIt.getIndices(i));
882 faceNodes[axisZ][faceIdxZ] = nS_Face;
884 newDataSourceIds->push_back(cellIt.getCorner(C).getPointId());
891 unsigned nNA = commitFace(axisA, A, P_A);
892 unsigned nNB = commitFace(axisB, B, P_B);
895 auto getEdgeNode = [&](
int v1,
int v2) {
897 for (
int e = 0; e < 12; ++e) {
898 if ((
corner0[e] ==
static_cast<unsigned>(v1) &&
899 corner1[e] ==
static_cast<unsigned>(v2)) ||
900 (
corner0[e] ==
static_cast<unsigned>(v2) &&
901 corner1[e] ==
static_cast<unsigned>(v1))) {
906 return this->
getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
909 unsigned nI_AD = getEdgeNode(A, D_corner);
910 unsigned nI_AZ = getEdgeNode(A, A_z);
911 unsigned nI_BD = getEdgeNode(B, D_corner);
912 unsigned nI_BZ = getEdgeNode(B, B_z);
913 unsigned nI_CZ = getEdgeNode(C, C_z);
916 int parity = (C & 1) + ((C >> 1) & 1) + ((C >> 2) & 1);
917 bool flip = (parity % 2 != 0) ^ inverted;
920 while ((A ^ C) != (1 << vA))
923 while ((B ^ C) != (1 << vB))
925 bool is_cyclic = ((vA + 1) % 3 == vB);
928 std::swap(nI_AD, nI_BD);
929 std::swap(nI_AZ, nI_BZ);
957 template <
int Dim = D>
959 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
960 int transform,
int axis,
bool inverted,
961 std::map<hrleIndex, unsigned> *nodes,
962 std::map<hrleIndex, unsigned> *faceNodes,
965 "Normal vector data must be computed for sharp edge generation.");
968 auto toCanonical = [&](Vec3D<T> v) {
984 res = Vec3D<T>{v[1], v[2], v[0]};
986 res = Vec3D<T>{v[2], v[0], v[1]};
991 auto fromCanonical = [&](Vec3D<T> v) {
997 res = Vec3D<T>{v[2], v[0], v[1]};
999 res = Vec3D<T>{v[1], v[2], v[0]};
1002 res[0] =
T(1.0) - res[0];
1004 res[1] =
T(1.0) - res[1];
1006 res[2] =
T(1.0) - res[2];
1010 auto getNormal = [&](
int idx) {
1011 auto corner = cellIt.getCorner(idx);
1012 if (corner.isDefined()) {
1013 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1015 for (
int i = 0; i < 3; ++i)
1025 unsigned faceNodeIds[2];
1026 bool nodeExists[2] = {
false,
false};
1029 auto mapIdx = [&](
int c) {
1032 int x = (c & 1), y = (c >> 1) & 1, z = (c >> 2) & 1;
1038 }
else if (axis == 1) {
1056 return gx | (gy << 1) | (gz << 2);
1060 for (
int k = 0; k < 2; ++k) {
1067 hrleIndex faceIdx = cellIt.getIndices();
1071 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1075 auto it = faceNodes[axis].find(faceIdx);
1076 if (it != faceNodes[axis].end()) {
1077 nodeExists[k] =
true;
1078 faceNodeIds[k] = it->second;
1079 Vec3D<T> relPos =
mesh->nodes[faceNodeIds[k]];
1080 for (
int d = 0; d < 3; ++d) {
1081 relPos[d] -=
static_cast<T>(cellIt.getIndices(d));
1083 P[k] = toCanonical(relPos);
1099 Vec3D<T> n_y = toCanonical(getNormal(c_y));
1100 Vec3D<T> n_z = toCanonical(getNormal(c_z));
1102 if (std::abs(DotProduct(n_y, n_z)) >= 0.8) {
1110 T t_y =
getInterp(c0, c_y, cellIt, inverted);
1111 T t_z =
getInterp(c0, c_z, cellIt, inverted);
1121 double det = n_y[1] * n_z[2] - n_y[2] * n_z[1];
1122 if (std::abs(det) < 1e-6)
1125 double d1 = n_y[1] * t_y;
1126 double d2 = n_z[2] * t_z;
1128 double y = (d1 * n_z[2] - d2 * n_y[2]) / det;
1129 double z = (d2 * n_y[1] - d1 * n_z[1]) / det;
1131 P[k] = Vec3D<T>{(
T)k, (
T)y, (
T)z};
1135 if (Norm2(P[k] - Vec3D<T>{(
T)k, t_y, 0}) > 9.0 ||
1136 Norm2(P[k] - Vec3D<T>{(
T)k, 0, t_z}) > 9.0) {
1143 for (
int k = 0; k < 2; ++k) {
1144 if (!nodeExists[k]) {
1146 hrleIndex faceIdx = cellIt.getIndices();
1147 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1151 Vec3D<T> globalP = fromCanonical(P[k]);
1152 for (
int d = 0; d < 3; ++d)
1153 globalP[d] = (globalP[d] + cellIt.getIndices(d));
1156 faceNodes[axis][faceIdx] = faceNodeIds[k];
1158 newDataSourceIds->push_back(
1159 cellIt.getCorner(mapIdx(c0)).getPointId());
1174 auto getEdgeNode = [&](
int c1,
int c2) {
1176 int g1 = mapIdx(c1);
1177 int g2 = mapIdx(c2);
1181 for (
int e = 0; e < 12; ++e) {
1182 if ((
corner0[e] ==
static_cast<unsigned>(g1) &&
1183 corner1[e] ==
static_cast<unsigned>(g2)) ||
1184 (
corner0[e] ==
static_cast<unsigned>(g2) &&
1185 corner1[e] ==
static_cast<unsigned>(g1))) {
1190 return this->
getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1193 unsigned n02 = getEdgeNode(0, 2);
1194 unsigned n04 = getEdgeNode(0, 4);
1195 unsigned n13 = getEdgeNode(1, 3);
1196 unsigned n15 = getEdgeNode(1, 5);
1214 template <
int Dim = D>
1216 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
int countNeg,
1217 int countPos,
int negMask,
int posMask,
1218 std::map<hrleIndex, unsigned> *nodes,
1219 std::map<hrleIndex, unsigned> *faceNodes,
1222 bool inverted =
false;
1224 if (countNeg == 2) {
1227 }
else if (countPos == 2) {
1235 int v1 = -1, v2 = -1;
1236 for (
int i = 0; i < 8; ++i) {
1237 if ((mask >> i) & 1) {
1246 if ((diff & (diff - 1)) != 0)
1250 while ((diff >> (axis + 1)) > 0)
1257 nodes, faceNodes, newDataSourceIds,
1264 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
1265 int transform,
bool inverted, std::map<hrleIndex, unsigned> *nodes,
1266 std::map<hrleIndex, unsigned> *faceNodes,
1268 Vec3D<T> &cornerPos) {
1274 "Normal vector data must be available for sharp corner generation");
1277 auto toCanonical = [&](Vec3D<T> v) {
1288 auto fromCanonical = [&](Vec3D<T> v) {
1290 v[0] =
T(1.0) - v[0];
1292 v[1] =
T(1.0) - v[1];
1294 v[2] =
T(1.0) - v[2];
1298 auto getNormal = [&](
int idx) {
1299 auto corner = cellIt.getCorner(idx);
1300 if (corner.isDefined()) {
1301 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1303 for (
int i = 0; i < 3; ++i)
1311 Vec3D<T> norm1, norm2, norm3;
1313 Vec3D<T> P1, P2, P4;
1315 const int g1 = transform ^ 1;
1316 const int g2 = transform ^ 2;
1317 const int g4 = transform ^ 4;
1319 norm1 = toCanonical(getNormal(g1));
1320 norm2 = toCanonical(getNormal(g2));
1321 norm3 = toCanonical(getNormal(g4));
1323 const T t1 =
getInterp(transform, g1, cellIt, inverted);
1324 const T t2 =
getInterp(transform, g2, cellIt, inverted);
1325 const T t4 =
getInterp(transform, g4, cellIt, inverted);
1331 d1 = DotProduct(norm1, P1);
1332 d2 = DotProduct(norm2, P2);
1333 d3 = DotProduct(norm3, P4);
1335 if (std::abs(DotProduct(norm1, norm2)) >= 0.8 ||
1336 std::abs(DotProduct(norm1, norm3)) >= 0.8 ||
1337 std::abs(DotProduct(norm2, norm3)) >= 0.8) {
1350 const Vec3D<T> c23 = CrossProduct(norm2, norm3);
1351 const Vec3D<T> c31 = CrossProduct(norm3, norm1);
1352 const Vec3D<T> c12 = CrossProduct(norm1, norm2);
1354 const double det = DotProduct(norm1, c23);
1359 const T invDet = 1.0 / det;
1361 for (
int i = 0; i < 3; ++i) {
1362 S[i] = (d1 * c23[i] + d2 * c31[i] + d3 * c12[i]) * invDet;
1365 if (Norm2(S - P1) > 9.0 || Norm2(S - P2) > 9.0 || Norm2(S - P4) > 9.0) {
1370 Vec3D<T> globalS = fromCanonical(S);
1371 for (
int d = 0; d < 3; ++d) {
1372 cornerPos[d] = (globalS[d] + cellIt.getIndices(d));
1380 template <
int Dim = D>
1382 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
int countNeg,
1383 int countPos,
int negMask,
int posMask,
1384 std::map<hrleIndex, unsigned> *nodes,
1385 std::map<hrleIndex, unsigned> &cornerNodes,
1386 std::map<hrleIndex, unsigned> *faceNodes,
1389 bool inverted =
false;
1393 if (countNeg == 1) {
1395 while (!((negMask >> transform) & 1))
1398 }
else if (countPos == 1) {
1400 while (!((posMask >> transform) & 1))
1405 if (transform == -1)
1410 faceNodes, newDataSourceIds, valueIt,
1413 hrleIndex cornerIdx = cellIt.getIndices();
1414 cornerIdx += viennahrle::BitMaskToIndex<D>(transform);
1416 auto it = cornerNodes.find(cornerIdx);
1417 if (it != cornerNodes.end()) {
1418 nCorner = it->second;
1421 cornerNodes[cornerIdx] = nCorner;
1427 newDataSourceIds->push_back(cellIt.getCorner(transform).getPointId());
1431 auto getEdgeNode = [&](
int neighbor) {
1435 for (
int e = 0; e < 12; ++e) {
1436 if ((
corner0[e] ==
static_cast<unsigned>(v1) &&
1437 corner1[e] ==
static_cast<unsigned>(v2)) ||
1438 (
corner0[e] ==
static_cast<unsigned>(v2) &&
1439 corner1[e] ==
static_cast<unsigned>(v1))) {
1444 return this->
getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1447 unsigned n1 = getEdgeNode(transform ^ 1);
1448 unsigned n2 = getEdgeNode(transform ^ 2);
1449 unsigned n4 = getEdgeNode(transform ^ 4);
1454 auto toCanonical = [&](Vec3D<T> v) {
1464 auto fromCanonical = [&](Vec3D<T> v) {
1466 v[0] =
T(1.0) - v[0];
1468 v[1] =
T(1.0) - v[1];
1470 v[2] =
T(1.0) - v[2];
1474 auto getNormal = [&](
int idx) {
1475 auto corner = cellIt.getCorner(idx);
1476 if (corner.isDefined()) {
1477 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1479 for (
int i = 0; i < 3; ++i)
1486 Vec3D<T> norm1, norm2, norm3;
1490 int n1_idx = -1, n2_idx = -1, n3_idx = -1;
1491 int neighbors[] = {transform ^ 1, transform ^ 2, transform ^ 4};
1493 for (
int n_idx : neighbors) {
1494 T val = cellIt.getCorner(n_idx).getValue();
1496 if (val * cellIt.getCorner(transform).getValue() < 0) {
1507 norm1 = toCanonical(getNormal(transform ^ 1));
1508 norm2 = toCanonical(getNormal(transform ^ 2));
1509 norm3 = toCanonical(getNormal(transform ^ 4));
1511 d1 = norm1[0] *
getInterp(transform, transform ^ 1, cellIt, inverted);
1512 d2 = norm2[1] *
getInterp(transform, transform ^ 2, cellIt, inverted);
1513 d3 = norm3[2] *
getInterp(transform, transform ^ 4, cellIt, inverted);
1515 auto solve2x2 = [&](
double a1,
double b1,
double c1,
double a2,
1516 double b2,
double c2,
T &res1,
T &res2) {
1517 double det = a1 * b2 - a2 * b1;
1518 if (std::abs(det) < 1e-6)
1520 res1 = (c1 * b2 - c2 * b1) / det;
1521 res2 = (a1 * c2 - a2 * c1) / det;
1525 Vec3D<T> Fx, Fy, Fz;
1531 valid &= solve2x2(norm2[1], norm2[2], d2 - norm2[0] * Fx[0], norm3[1],
1532 norm3[2], d3 - norm3[0] * Fx[0], Fx[1], Fx[2]);
1533 valid &= solve2x2(norm1[0], norm1[2], d1 - norm1[1] * Fy[1], norm3[0],
1534 norm3[2], d3 - norm3[1] * Fy[1], Fy[0], Fy[2]);
1535 valid &= solve2x2(norm1[0], norm1[1], d1 - norm1[2] * Fz[2], norm2[0],
1536 norm2[1], d2 - norm2[2] * Fz[2], Fz[0], Fz[1]);
1538 auto checkBounds = [&](Vec3D<T> p) {
1539 return p[0] >= -1e-4 && p[0] <= 1.0 + 1e-4 && p[1] >= -1e-4 &&
1540 p[1] <= 1.0 + 1e-4 && p[2] >= -1e-4 && p[2] <= 1.0 + 1e-4;
1543 if (valid && checkBounds(Fx) && checkBounds(Fy) && checkBounds(Fz)) {
1544 auto getFaceNode = [&](
int axis, Vec3D<T> pos) {
1545 hrleIndex faceIdx = cellIt.getIndices();
1546 if ((transform >> axis) & 1)
1549 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
1550 return faceNodes[axis][faceIdx];
1552 Vec3D<T> globalPos = fromCanonical(pos);
1553 for (
int d = 0; d < 3; ++d)
1554 globalPos[d] = (globalPos[d] + cellIt.getIndices(d));
1556 faceNodes[axis][faceIdx] = nodeId;
1558 newDataSourceIds->push_back(
1559 cellIt.getCorner(transform).getPointId());
1565 unsigned nFx = getFaceNode(0, Fx);
1566 unsigned nFy = getFaceNode(1, Fy);
1567 unsigned nFz = getFaceNode(2, Fz);
1571 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1574 bool flip = (parity % 2 != 0) ^ inverted;
1602 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1605 bool flip = (parity % 2 != 0) ^ inverted;
1625 if constexpr (
D == 3) {
1626 return nodeNumbers[0] == nodeNumbers[1] ||
1627 nodeNumbers[0] == nodeNumbers[2] ||
1628 nodeNumbers[1] == nodeNumbers[2];
1630 return nodeNumbers[0] == nodeNumbers[1];
1635 const Vec3D<T> &nodeB,
1636 const Vec3D<T> &nodeC)
noexcept {
1637 return CrossProduct(nodeB - nodeA, nodeC - nodeA);
1641 const viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
1642 bool inverted)
const {
1643 auto getValue = [&](
int idx) {
1644 T val = cellIt.getCorner(idx).getValue();
1645 return inverted ? -val : val;
1647 const T v_a = getValue(p_a);
1648 const T v_b = getValue(p_b);
1649 if (std::abs(v_a - v_b) <
epsilon)
1651 return v_a / (v_a - v_b);
1658 auto elementI3 = [&](
const std::array<unsigned, D> &element) ->
I3 {
1659 if constexpr (
D == 2)
1660 return {(int)element[0], (
int)element[1], 0};
1662 return {(int)element[0], (
int)element[1], (int)element[2]};
1667 if constexpr (
D == 2) {
1669 -(
mesh->nodes[nodeNumbers[1]][1] -
mesh->nodes[nodeNumbers[0]][1]),
1670 mesh->nodes[nodeNumbers[1]][0] -
mesh->nodes[nodeNumbers[0]][0],
1674 mesh->nodes[nodeNumbers[1]],
1675 mesh->nodes[nodeNumbers[2]]);
1679 normal[0] * normal[0] + normal[1] * normal[1] + normal[2] * normal[2];
1681 mesh->insertNextElement(nodeNumbers);
1682 T invn =
static_cast<T>(1.) / std::sqrt(
static_cast<T>(n2));
1683 for (
int d = 0; d <
D; d++) {
1693 auto quantize = [&](
const Vec3D<T> &p) ->
I3 {
1695 return {(int)std::llround(p[0] * inv), (int)std::llround(p[1] * inv),
1696 (int)std::llround(p[2] * inv)};
1701 auto q = quantize(pos);
1704 nodeIdx = it->second;
1708 for (
int i = 0; i < 3; ++i) {
1709 mesh->nodes[nodeIdx][i] = (
mesh->nodes[nodeIdx][i] + pos[i]) *
T(0.5);
1713 unsigned newNodeId =
mesh->insertNextNode(pos);
1717 for (
int a = 0; a <
D; a++) {
1718 if (pos[a] <
mesh->minimumExtent[a])
1719 mesh->minimumExtent[a] = pos[a];
1720 if (pos[a] >
mesh->maximumExtent[a])
1721 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
std::enable_if_t< Dim==3, 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:675
std::enable_if_t< Dim==3, 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:958
static constexpr unsigned int direction[12]
Definition lsToSurfaceMesh.hpp:77
void setMesh(SmartPointer< Mesh< T > > passedMesh)
Definition lsToSurfaceMesh.hpp:95
SmartPointer< lsDomainType > currentLevelSet
Definition lsToSurfaceMesh.hpp:67
std::enable_if_t< Dim==3, 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:1215
bool generateSharpCorners
Definition lsToSurfaceMesh.hpp:38
std::vector< SmartPointer< lsDomainType > > levelSets
Definition lsToSurfaceMesh.hpp:32
std::enable_if_t< Dim==2, bool > generateCanonicalSharpCorner2D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int transform, Vec3D< T > &cornerPos) const
Definition lsToSurfaceMesh.hpp:480
typename lsDomainType::DomainType hrleDomainType
Definition lsToSurfaceMesh.hpp:28
static bool triangleMisformed(const std::array< unsigned, D > &nodeNumbers) noexcept
Definition lsToSurfaceMesh.hpp:1624
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
std::unordered_set< I3, I3Hash > uniqueElements
Definition lsToSurfaceMesh.hpp:62
std::enable_if_t< Dim==3, 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:1381
std::enable_if_t< Dim==3, void > stitchToNeighbor(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int axis, bool isHighFace, unsigned faceNodeId, std::map< hrleIndex, unsigned > *nodes, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:423
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:1634
const double minNodeDistanceFactor
Definition lsToSurfaceMesh.hpp:36
void insertElement(const std::array< unsigned, D > &nodeNumbers)
Definition lsToSurfaceMesh.hpp:1654
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
std::enable_if_t< Dim==2, 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:555
unsigned insertNode(Vec3D< T > const &pos)
Definition lsToSurfaceMesh.hpp:1692
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:1263
T getInterp(int p_a, int p_b, const viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, bool inverted) const
Definition lsToSurfaceMesh.hpp:1640
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