ViennaLS
Loading...
Searching...
No Matches
lsToSurfaceMesh.hpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <map>
6#include <unordered_map>
7#include <unordered_set>
8
9#include <hrleSparseCellIterator.hpp>
10#include <hrleSparseStarIterator.hpp>
12#include <lsDomain.hpp>
13#include <lsExpand.hpp>
15#include <lsMarchingCubes.hpp>
16#include <lsMesh.hpp>
17
18namespace viennals {
19
20using namespace viennacore;
21
25template <class T, int D> class ToSurfaceMesh {
26protected:
29 using hrleIndex = viennahrle::Index<D>;
30 using ConstSparseIterator = viennahrle::ConstSparseIterator<hrleDomainType>;
31
32 std::vector<SmartPointer<lsDomainType>> levelSets;
33 SmartPointer<Mesh<T>> mesh = nullptr;
34
35 const T epsilon;
37 bool updatePointData = true;
39
40 struct I3 {
41 int x, y, z;
42 bool operator==(const I3 &o) const {
43 return x == o.x && y == o.y && z == o.z;
44 }
45 };
46
47 struct I3Hash {
48 size_t operator()(const I3 &k) const {
49 // 64-bit mix
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);
56 return (size_t)h;
57 }
58 };
59
60 // Context for mesh generation to share between helpers
61 std::unordered_map<I3, unsigned, I3Hash> nodeIdByBin;
62 std::unordered_set<I3, I3Hash> uniqueElements;
63 std::vector<Vec3D<T>> currentNormals;
64 std::vector<T> currentMaterials;
67 SmartPointer<lsDomainType> currentLevelSet = nullptr;
69
70 // Store sharp corner nodes created during this cell iteration
71 std::vector<std::pair<unsigned, Vec3D<T>>> matSharpCornerNodes;
72
73 static constexpr unsigned int corner0[12] = {0, 1, 2, 0, 4, 5,
74 6, 4, 0, 1, 3, 2};
75 static constexpr unsigned int corner1[12] = {1, 3, 3, 2, 5, 7,
76 7, 6, 4, 5, 7, 6};
77 static constexpr unsigned int direction[12] = {0, 1, 0, 1, 0, 1,
78 0, 1, 2, 2, 2, 2};
79
80public:
81 explicit ToSurfaceMesh(double mnd = 0.05, double eps = 1e-12)
82 : minNodeDistanceFactor(mnd), epsilon(eps) {}
83
84 ToSurfaceMesh(const SmartPointer<lsDomainType> passedLevelSet,
85 SmartPointer<Mesh<T>> passedMesh, double mnd = 0.05,
86 double eps = 1e-12)
87 : levelSets({passedLevelSet}), mesh(passedMesh),
88 minNodeDistanceFactor(mnd), epsilon(eps) {}
89
90 void setLevelSet(SmartPointer<lsDomainType> passedlsDomain) {
91 levelSets = {passedlsDomain};
92 normalVectorData = nullptr;
93 }
94
95 void setMesh(SmartPointer<Mesh<T>> passedMesh) { mesh = passedMesh; }
96
97 void setUpdatePointData(bool update) { updatePointData = update; }
98
99 void setSharpCorners(bool check) { generateSharpCorners = check; }
100
101 virtual void apply() {
103 currentGridDelta = currentLevelSet->getGrid().getGridDelta();
104 if (currentLevelSet == nullptr) {
105 Logger::getInstance()
106 .addError("No level set was passed to ToSurfaceMesh.")
107 .print();
108 return;
109 }
110 if (mesh == nullptr) {
111 Logger::getInstance()
112 .addError("No mesh was passed to ToSurfaceMesh.")
113 .print();
114 return;
115 }
116
117 if (currentLevelSet->getNumberOfPoints() == 0) {
118 VIENNACORE_LOG_WARNING(
119 "ToSurfaceMesh: Level set is empty. No mesh will be created.");
120 return;
121 }
122
123 mesh->clear();
124 if (currentLevelSet != nullptr) {
125 if (currentLevelSet->getLevelSetWidth() < 2) {
126 VIENNACORE_LOG_WARNING("Levelset is less than 2 layers wide. Expanding "
127 "levelset to 2 layers.");
129 }
130 }
131
132 ConstSparseIterator valueIt(currentLevelSet->getDomain());
133
134 typedef std::map<hrleIndex, unsigned> nodeContainerType;
135
136 nodeContainerType nodes[D];
137 nodeContainerType faceNodes[D];
138 nodeContainerType cornerNodes;
139 const bool sharpCorners = generateSharpCorners;
140 typename nodeContainerType::iterator nodeIt;
141
142 // save how data should be transferred to new level set
143 // list of indices into the old pointData vector
144 std::vector<std::vector<unsigned>> newDataSourceIds;
145 // there is no multithreading here, so just use 1
146 if (updatePointData)
147 newDataSourceIds.resize(1);
148
149 // If sharp corners are enabled, calculate normals first as they are needed
150 // for feature reconstruction
151 if (sharpCorners) {
153 normalCalculator.setMethod(
155 normalCalculator.setMaxValue(std::numeric_limits<T>::max());
156 normalCalculator.apply();
157 normalVectorData = currentLevelSet->getPointData().getVectorData(
159 }
160
161 // iterate over all cells with active points
162 for (viennahrle::ConstSparseCellIterator<hrleDomainType> cellIt(
163 currentLevelSet->getDomain());
164 !cellIt.isFinished(); cellIt.next()) {
165
166 // Clear node caches for dimensions that have moved out of scope
167 if (!sharpCorners) {
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());
172
173 while (!faceNodes[u].empty() &&
174 faceNodes[u].begin()->first < hrleIndex(cellIt.getIndices()))
175 faceNodes[u].erase(faceNodes[u].begin());
176
177 if (u == 0) {
178 while (!cornerNodes.empty() &&
179 cornerNodes.begin()->first < hrleIndex(cellIt.getIndices()))
180 cornerNodes.erase(cornerNodes.begin());
181 }
182 }
183 }
184
185 // Calculate signs of all corners to determine the marching cubes case
186 unsigned signs = 0;
187 bool hasZero = false;
188 for (int i = 0; i < (1 << D); i++) {
189 T val = cellIt.getCorner(i).getValue();
190 if (val >= T(0))
191 signs |= (1 << i);
192 if (std::abs(val) <= epsilon)
193 hasZero = true;
194 }
195
196 // all corners have the same sign, so no surface here
197 if (signs == 0)
198 continue;
199 if (signs == (1 << (1 << D)) - 1 && !hasZero)
200 continue;
201
202 // Attempt to generate sharp features if enabled
203 bool perfectCornerFound = false;
204 if (sharpCorners) {
205 int countNeg = 0;
206 int countPos = 0;
207 int negMask = 0;
208 int posMask = 0;
209 for (int i = 0; i < (1 << D); ++i) {
210 T val = cellIt.getCorner(i).getValue();
211 if (val < -epsilon) {
212 countNeg++;
213 negMask |= (1 << i);
214 } else if (val > epsilon) {
215 countPos++;
216 posMask |= (1 << i);
217 } else {
218 if (val >= 0) {
219 countPos++;
220 posMask |= (1 << i);
221 } else {
222 countNeg++;
223 negMask |= (1 << i);
224 }
225 }
226 }
227
228 // Check for perfect corner (2D only)
229 if constexpr (D == 2) {
230 perfectCornerFound = generateSharpCorner2D(
231 cellIt, countNeg, countPos, negMask, posMask, nodes,
232 &newDataSourceIds[0], valueIt);
233 } else if constexpr (D == 3) {
234 if (countNeg == 2 || countPos == 2) {
235 // Try to generate a sharp edge (2 corners active)
236 perfectCornerFound = generateSharpEdge3D(
237 cellIt, countNeg, countPos, negMask, posMask, nodes, faceNodes,
238 &newDataSourceIds[0], valueIt);
239 } else if (countNeg == 1 || countPos == 1) {
240 // Try to generate a sharp corner (1 corner active)
241 perfectCornerFound = generateSharpCorner3D(
242 cellIt, countNeg, countPos, negMask, posMask, nodes,
243 cornerNodes, faceNodes, &newDataSourceIds[0], valueIt);
244 } else if (countNeg == 3 || countPos == 3) {
245 // Try to generate an L-shape corner (3 corners active)
246 perfectCornerFound = generateSharpL3D(
247 cellIt, countNeg, countPos, negMask, posMask, nodes, faceNodes,
248 &newDataSourceIds[0], valueIt);
249 }
250 }
251 }
252
253 // If a sharp feature was successfully generated, skip standard Marching
254 // Cubes for this cell
255 if (perfectCornerFound)
256 continue;
257
258 if constexpr (D == 3) {
259 // Stitch to perfect corners/edges
260 // Check if neighbors have generated nodes on shared faces that we need
261 // to connect to
262 for (int axis = 0; axis < 3; ++axis) {
263 for (int d = 0; d < 2; ++d) {
264 hrleIndex faceKey(cellIt.getIndices());
265 if (d == 1)
266 faceKey[axis]++;
267
268 auto it = faceNodes[axis].find(faceKey);
269 if (it != faceNodes[axis].end()) {
270 const unsigned faceNodeId = it->second;
271 const int *Triangles =
273
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];
278 int c0 = corner0[edge];
279 int c1 = corner1[edge];
280 bool onFace =
281 (((c0 >> axis) & 1) == d) && (((c1 >> axis) & 1) == d);
282 if (onFace) {
283 face_edge_nodes.push_back(
284 getNode(cellIt, edge, nodes, &newDataSourceIds[0]));
285 }
286 }
287 if (face_edge_nodes.size() == 2) {
289 {face_edge_nodes[0], face_edge_nodes[1], faceNodeId});
290 }
291 }
292 }
293 }
294 }
295 }
296
297 // Standard Marching Cubes / Squares algorithm
298 // for each element
299 const int *Triangles =
302
303 for (; Triangles[0] != -1; Triangles += D) {
304 std::array<unsigned, D> nodeNumbers;
305
306 // for each node
307 for (int n = 0; n < D; n++) {
308 const int edge = Triangles[n];
309
310 unsigned p0 = corner0[edge];
311 unsigned p1 = corner1[edge];
312
313 // determine direction of edge
314 auto dir = direction[edge];
315
316 // look for existing surface node
317 hrleIndex d(cellIt.getIndices());
318 d += viennahrle::BitMaskToIndex<D>(p0);
319
320 nodeIt = nodes[dir].find(d);
321 if (nodeIt != nodes[dir].end()) {
322 nodeNumbers[n] = nodeIt->second;
323 } else { // if node does not exist yet
324 nodeNumbers[n] = getNode(cellIt, edge, nodes, &newDataSourceIds[0]);
325 }
326 }
327
328 insertElement(nodeNumbers); // insert new surface element
329 }
330 }
331
332 scaleMesh();
333
334 // now copy old data into new level set
335 if (updatePointData && !newDataSourceIds[0].empty()) {
336 mesh->getPointData().translateFromMultiData(
337 currentLevelSet->getPointData(), newDataSourceIds);
338 }
339 }
340
341protected:
342 void scaleMesh() {
343 // Scale the mesh to global coordinates
344 for (auto &node : mesh->nodes) {
345 for (int i = 0; i < D; ++i)
346 node[i] *= currentGridDelta;
347 }
348 for (int i = 0; i < D; ++i) {
349 mesh->minimumExtent[i] *= currentGridDelta;
350 mesh->maximumExtent[i] *= currentGridDelta;
351 }
352 }
353
354 // Helper to get or create a node on an edge using linear interpolation
355 // Compute the interpolated position of a node on an edge without inserting it
357 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt, int edge) {
358 unsigned p0 = corner0[edge];
359 unsigned p1 = corner1[edge];
360 auto dir = direction[edge];
361
362 Vec3D<T> cc{};
363 for (int z = 0; z < D; z++) {
364 if (z != dir) {
365 cc[z] = static_cast<double>(cellIt.getIndices(z) +
366 viennahrle::BitMaskToIndex<D>(p0)[z]);
367 } else {
368 T d0 = cellIt.getCorner(p0).getValue();
369 T d1 = cellIt.getCorner(p1).getValue();
370 if (d0 == -d1) {
371 cc[z] = static_cast<T>(cellIt.getIndices(z)) + T(0.5);
372 } else {
373 if (std::abs(d0) <= std::abs(d1)) {
374 cc[z] = static_cast<T>(cellIt.getIndices(z)) + (d0 / (d0 - d1));
375 } else {
376 cc[z] = static_cast<T>(cellIt.getIndices(z) + 1) - (d1 / (d1 - d0));
377 }
378 }
379 cc[z] = std::max(cc[z], cellIt.getIndices(z) + epsilon);
380 cc[z] = std::min(cc[z], (cellIt.getIndices(z) + 1) - epsilon);
381 }
382 }
383 return cc;
384 }
385
386 unsigned getNode(viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
387 int edge, std::map<hrleIndex, unsigned> *nodes,
388 std::vector<unsigned> *newDataSourceIds) {
389
390 unsigned p0 = corner0[edge];
391 unsigned p1 = corner1[edge];
392 auto dir = direction[edge];
393 hrleIndex d(cellIt.getIndices());
394 d += viennahrle::BitMaskToIndex<D>(p0);
395
396 if (nodes[dir].find(d) != nodes[dir].end()) {
397 return nodes[dir][d];
398 }
399
400 // Create node
401 Vec3D<T> cc = computeNodePosition(cellIt, edge);
402
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();
408 } else {
409 currentPointId = cellIt.getCorner(p1).getPointId();
410 }
411 if (updatePointData && newDataSourceIds)
412 newDataSourceIds->push_back(currentPointId);
413
414 unsigned nodeId = insertNode(cc);
415 nodes[dir][d] = nodeId;
416 return nodeId;
417 }
418
419 // Helper to stitch a sharp feature node on a face to the standard mesh in a
420 // neighboring cell (3D only)
421 template <int Dim = D>
422 std::enable_if_t<Dim == 3, void>
423 stitchToNeighbor(viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
424 int axis, bool isHighFace, unsigned faceNodeId,
425 std::map<hrleIndex, unsigned> *nodes,
426 ConstSparseIterator &valueIt) {
427 // Backward stitching: Check if neighbor on this face is "past" and needs
428 // stitching
429 hrleIndex neighborIdx = cellIt.getIndices();
430 if (isHighFace)
431 neighborIdx[axis]++;
432 else
433 neighborIdx[axis]--;
434
435 if (neighborIdx < cellIt.getIndices()) {
436 unsigned nSigns = 0;
437 auto &grid = currentLevelSet->getGrid();
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)
443 nSigns |= (1 << i);
444 }
445 }
446
447 if (nSigns != 0 && nSigns != 255) {
448 const int *nTriangles = lsInternal::MarchingCubes::polygonize3d(nSigns);
449 int nFaceD = isHighFace ? 0 : 1;
450
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];
455 int c0 = corner0[edge];
456 int c1 = corner1[edge];
457 bool onFace = (((c0 >> axis) & 1) == nFaceD) &&
458 (((c1 >> axis) & 1) == nFaceD);
459 if (onFace) {
460 unsigned p0 = corner0[edge];
461 auto dir = direction[edge];
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);
466 }
467 }
468 }
469 if (face_edge_nodes.size() == 2) {
470 insertElement(std::array<unsigned, 3>{
471 face_edge_nodes[0], face_edge_nodes[1], faceNodeId});
472 }
473 }
474 }
475 }
476 }
477
478 // Solves for a sharp corner position in 2D by intersecting normals
479 template <int Dim = D>
480 std::enable_if_t<Dim == 2, bool> generateCanonicalSharpCorner2D(
481 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
482 int transform, Vec3D<T> &cornerPos) const {
483 assert(normalVectorData &&
484 "Normal vector data must be computed for sharp corner generation.");
485
486 auto getTransformedGradient = [&](int cornerID) {
487 auto corner = cellIt.getCorner(cornerID ^ transform);
488 if (corner.isDefined()) {
489 auto normal = (*normalVectorData)[corner.getPointId()];
490 normal[2] = 0;
491 if ((transform & 1) != 0)
492 normal[0] = -normal[0];
493 if ((transform & 2) != 0)
494 normal[1] = -normal[1];
495 return normal;
496 }
497 return Vec3D<T>{};
498 };
499
500 Vec3D<T> norm1 = getTransformedGradient(1); // neighbor in x
501 Vec3D<T> norm2 = getTransformedGradient(2); // neighbor in y
502
503 if (std::abs(DotProduct(norm1, norm2)) >= 0.8) {
504 return false;
505 }
506
507 auto calculateNodePos = [&](int edge, Vec3D<T> &pos) {
508 unsigned p0 = corner0[edge];
509 unsigned p1 = corner1[edge];
510 auto dir = direction[edge];
511
512 for (int z = 0; z < D; z++) {
513 if (z != dir) {
514 pos[z] = static_cast<T>(viennahrle::BitMaskToIndex<D>(p0)[z]);
515 } else {
516 const T d0 = cellIt.getCorner(p0 ^ transform).getValue();
517 const T d1 = cellIt.getCorner(p1 ^ transform).getValue();
518 if (d0 == -d1) {
519 pos[z] = T(0.5);
520 } else {
521 if (std::abs(d0) <= std::abs(d1)) {
522 pos[z] = (d0 / (d0 - d1));
523 } else {
524 pos[z] = T(1.0) - (d1 / (d1 - d0));
525 }
526 }
527 pos[z] = std::max(pos[z], epsilon);
528 pos[z] = std::min(pos[z], T(1.0) - epsilon);
529 }
530 }
531 };
532
533 Vec3D<T> pX{}, pY{};
534 calculateNodePos(0, pX); // Edge along x-axis from corner 0
535 calculateNodePos(3, pY); // Edge along y-axis from corner 0
536
537 double d1 = DotProduct(norm1, pX);
538 double d2 = DotProduct(norm2, pY);
539 double det = norm1[0] * norm2[1] - norm1[1] * norm2[0];
540
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;
544 } else {
545 for (int i = 0; i < D; ++i)
546 cornerPos[i] = (pX[i] + pY[i]) * T(0.5);
547 }
548
549 return true;
550 }
551
552 // Wrapper for 2D sharp corner generation handling different
553 // rotations/reflections
554 template <int Dim = D>
555 std::enable_if_t<Dim == 2, bool> generateSharpCorner2D(
556 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt, int countNeg,
557 int countPos, int negMask, int posMask,
558 std::map<hrleIndex, unsigned> *nodes,
559 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
560
561 int cornerIdx = -1;
562 if (countNeg == 1) {
563 for (int i = 0; i < 4; ++i)
564 if ((negMask >> i) & 1) {
565 cornerIdx = i;
566 break;
567 }
568 } else if (countPos == 1) {
569 for (int i = 0; i < 4; ++i)
570 if ((posMask >> i) & 1) {
571 cornerIdx = i;
572 break;
573 }
574 }
575
576 if (cornerIdx != -1) {
577 // Check if this corner is also a corner in any neighboring cell
578 bool isSharedCorner = false;
579 auto pIdx =
580 cellIt.getIndices() + viennahrle::BitMaskToIndex<D>(cornerIdx);
581 T pVal = cellIt.getCorner(cornerIdx).getValue();
582 auto &grid = currentLevelSet->getGrid();
583
584 for (int i = 0; i < (1 << D); ++i) {
585 if (i == cornerIdx)
586 continue;
587
588 hrleIndex neighborIndices;
589 for (int k = 0; k < D; ++k)
590 neighborIndices[k] = pIdx[k] - ((i >> k) & 1);
591
592 bool neighborIsCorner = true;
593
594 // Check edge-connected neighbors in the neighbor cell
595 for (int k = 0; k < D; ++k) {
596 int neighborLocal = i ^ (1 << k);
597 auto checkIdx =
598 neighborIndices + viennahrle::BitMaskToIndex<D>(neighborLocal);
599 if (grid.isOutsideOfDomain(checkIdx)) {
600 checkIdx = grid.globalIndices2LocalIndices(checkIdx);
601 }
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;
608 break;
609 }
610 }
611
612 if (neighborIsCorner) {
613 isSharedCorner = true;
614 break;
615 }
616 }
617
618 if (!isSharedCorner) {
619 Vec3D<T> cornerPos{};
620 if (generateCanonicalSharpCorner2D(cellIt, cornerIdx, cornerPos)) {
621
622 // inverse transform cornerPos from canonical local to this cell's
623 // local
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];
628
629 // convert to global coordinates
630 for (int i = 0; i < D; ++i) {
631 cornerPos[i] = (cornerPos[i] + cellIt.getIndices(i));
632 }
633
634 // Determine edges for this corner
635 int edgeX = -1, edgeY = -1;
636 if (cornerIdx == 0) {
637 edgeX = 0;
638 edgeY = 3;
639 } else if (cornerIdx == 1) {
640 edgeX = 0;
641 edgeY = 1;
642 } else if (cornerIdx == 2) {
643 edgeX = 2;
644 edgeY = 3;
645 } else if (cornerIdx == 3) {
646 edgeX = 2;
647 edgeY = 1;
648 }
649
650 unsigned nX = getNode(cellIt, edgeX, nodes, newDataSourceIds);
651 unsigned nY = getNode(cellIt, edgeY, nodes, newDataSourceIds);
652 unsigned nCorner = insertNode(cornerPos);
653
654 // Store this sharp corner node for potential use by derived classes
655 matSharpCornerNodes.push_back({nCorner, cornerPos});
656
657 if (updatePointData && newDataSourceIds) {
658 assert(!newDataSourceIds->empty());
659 newDataSourceIds->push_back(
660 newDataSourceIds->back()); // TODO: improve point data source
661 }
662
663 insertElement({nX, nCorner});
664 insertElement({nCorner, nY});
665 return true;
666 }
667 }
668 }
669 return false;
670 }
671
672 // Generates geometry for an "L-shaped" configuration (3 active corners) in 3D
673 template <int Dim = D>
674 std::enable_if_t<Dim == 3, bool>
675 generateSharpL3D(viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
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,
680 ConstSparseIterator &valueIt) {
681
682 bool inverted = false;
683 int mask = 0;
684 if (countNeg == 3) {
685 mask = negMask;
686 inverted = false;
687 } else if (countPos == 3) {
688 mask = posMask;
689 inverted = true;
690 } else {
691 return false;
692 }
693
694 // Check for L-shape: 3 points on a face
695 int C = -1, A = -1, B = -1;
696 for (int i = 0; i < 8; ++i) {
697 if ((mask >> i) & 1) {
698 int neighbors = 0;
699 for (int k = 0; k < 3; ++k) {
700 if ((mask >> (i ^ (1 << k))) & 1)
701 neighbors++;
702 }
703 if (neighbors == 2) {
704 C = i;
705 break;
706 }
707 }
708 }
709
710 if (C == -1)
711 return false;
712
713 // Identify A and B
714 for (int k = 0; k < 3; ++k) {
715 int n = C ^ (1 << k);
716 if ((mask >> n) & 1) {
717 if (A == -1)
718 A = n;
719 else
720 B = n;
721 }
722 }
723
724 assert(normalVectorData &&
725 "Normal vector data must be computed for sharp feature generation.");
726
727 auto getNormal = [&](int idx) {
728 auto corner = cellIt.getCorner(idx);
729 if (corner.isDefined()) {
730 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
731 if (inverted)
732 for (int i = 0; i < 3; ++i)
733 n[i] = -n[i];
734 return n;
735 }
736 return Vec3D<T>{};
737 };
738
739 // Determine axes
740 int axisA = 0;
741 while ((C ^ A) != (1 << axisA))
742 axisA++;
743 int axisB = 0;
744 while ((C ^ B) != (1 << axisB))
745 axisB++;
746 int axisZ = 3 - axisA - axisB;
747
748 // Solve 2D face problem to find sharp point on a face
749 auto calculateFace = [&](int axis, int corner, int n1,
750 int n2) -> std::optional<Vec3D<T>> {
751 hrleIndex faceIdx = cellIt.getIndices();
752 if ((corner >> axis) & 1)
753 faceIdx[axis]++;
754
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);
760 return pos;
761 }
762
763 Vec3D<T> N1 = getNormal(n1);
764 Vec3D<T> N2 = getNormal(n2);
765
766 if (std::abs(DotProduct(N1, N2)) >= 0.8)
767 return std::nullopt;
768
769 int u = (axis + 1) % 3;
770 int v = (axis + 2) % 3;
771
772 int axis1 = 0;
773 while ((corner ^ n1) != (1 << axis1))
774 axis1++;
775 T t1 = getInterp(corner, n1, cellIt, inverted);
776
777 int axis2 = 0;
778 while ((corner ^ n2) != (1 << axis2))
779 axis2++;
780 T t2 = getInterp(corner, n2, cellIt, inverted);
781
782 Vec3D<T> P1;
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;
787
788 Vec3D<T> P2;
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;
793
794 double det = N1[u] * N2[v] - N1[v] * N2[u];
795 if (std::abs(det) < 1e-6)
796 return std::nullopt;
797
798 double c1 = N1[u] * P1[u] + N1[v] * P1[v];
799 double c2 = N2[u] * P2[u] + N2[v] * P2[v];
800
801 Vec3D<T> P;
802 P[axis] = P1[axis];
803 P[u] = (c1 * N2[v] - c2 * N1[v]) / det;
804 P[v] = (c2 * N1[u] - c1 * N2[u]) / det;
805
806 if (Norm2(P - P1) > 9.0 || Norm2(P - P2) > 9.0)
807 return std::nullopt;
808
809 return P;
810 };
811
812 auto commitFace = [&](int axis, int corner, Vec3D<T> P) -> unsigned {
813 hrleIndex faceIdx = cellIt.getIndices();
814 if ((corner >> axis) & 1)
815 faceIdx[axis]++;
816
817 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
818 return faceNodes[axis][faceIdx];
819 }
820
821 Vec3D<T> globalP = P;
822 for (int d = 0; d < 3; ++d)
823 globalP[d] = (globalP[d] + cellIt.getIndices(d));
824 unsigned nodeId = insertNode(globalP);
825 faceNodes[axis][faceIdx] = nodeId;
826 if (updatePointData && newDataSourceIds)
827 newDataSourceIds->push_back(cellIt.getCorner(corner).getPointId());
828
829 stitchToNeighbor(cellIt, axis, (corner >> axis) & 1, nodeId, nodes,
830 valueIt);
831
832 return nodeId;
833 };
834
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);
839
840 auto P_A_opt = calculateFace(axisA, A, D_corner, A_z);
841 auto P_B_opt = calculateFace(axisB, B, D_corner, B_z);
842
843 if (!P_A_opt.has_value() || !P_B_opt.has_value())
844 return false;
845
846 auto const &P_A = P_A_opt.value();
847 auto const &P_B = P_B_opt.value();
848
849 // Construct S
850 Vec3D<T> S;
851 S[axisA] = P_B[axisA];
852 S[axisB] = P_A[axisB];
853 S[axisZ] = (P_A[axisZ] + P_B[axisZ]) * 0.5;
854
855 unsigned nS = insertNode([&]() {
856 Vec3D<T> gS = S;
857 for (int i = 0; i < 3; ++i)
858 gS[i] = (gS[i] + cellIt.getIndices(i));
859 return gS;
860 }());
861 if (updatePointData && newDataSourceIds)
862 newDataSourceIds->push_back(cellIt.getCorner(C).getPointId());
863
864 // Calculate S_Face
865 Vec3D<T> S_Face = S;
866 S_Face[axisZ] = static_cast<T>((C >> axisZ) & 1);
867
868 unsigned nS_Face;
869 auto faceIdxZ = cellIt.getIndices();
870 if ((C >> axisZ) & 1)
871 faceIdxZ[axisZ]++;
872
873 if (faceNodes[axisZ].find(faceIdxZ) != faceNodes[axisZ].end()) {
874 nS_Face = faceNodes[axisZ][faceIdxZ];
875 } else {
876 nS_Face = insertNode([&]() {
877 Vec3D<T> gS = S_Face;
878 for (int i = 0; i < 3; ++i)
879 gS[i] = (gS[i] + cellIt.getIndices(i));
880 return gS;
881 }());
882 faceNodes[axisZ][faceIdxZ] = nS_Face;
883 if (updatePointData && newDataSourceIds)
884 newDataSourceIds->push_back(cellIt.getCorner(C).getPointId());
885
886 stitchToNeighbor(cellIt, axisZ, (C >> axisZ) & 1, nS_Face, nodes,
887 valueIt);
888 }
889
890 // Get face nodes
891 unsigned nNA = commitFace(axisA, A, P_A);
892 unsigned nNB = commitFace(axisB, B, P_B);
893
894 // Get boundary intersection nodes
895 auto getEdgeNode = [&](int v1, int v2) {
896 int edgeIdx = -1;
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))) {
902 edgeIdx = e;
903 break;
904 }
905 }
906 return this->getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
907 };
908
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);
914
915 // Winding order
916 int parity = (C & 1) + ((C >> 1) & 1) + ((C >> 2) & 1);
917 bool flip = (parity % 2 != 0) ^ inverted;
918
919 int vA = 0;
920 while ((A ^ C) != (1 << vA))
921 vA++;
922 int vB = 0;
923 while ((B ^ C) != (1 << vB))
924 vB++;
925 bool is_cyclic = ((vA + 1) % 3 == vB);
926 if (!is_cyclic) {
927 std::swap(nNA, nNB);
928 std::swap(nI_AD, nI_BD);
929 std::swap(nI_AZ, nI_BZ);
930 }
931
932 if (!flip) {
933 insertElement({nS, nNA, nI_AD});
934 insertElement({nS, nI_AD, nS_Face});
935 insertElement({nS, nS_Face, nI_BD});
936 insertElement({nS, nI_BD, nNB});
937 insertElement({nS, nNB, nI_BZ});
938 insertElement({nS, nI_BZ, nI_CZ});
939 insertElement({nS, nI_CZ, nI_AZ});
940 insertElement({nS, nI_AZ, nNA});
941 } else {
942 insertElement({nS, nI_AD, nNA});
943 insertElement({nS, nS_Face, nI_AD});
944 insertElement({nS, nI_BD, nS_Face});
945 insertElement({nS, nNB, nI_BD});
946 insertElement({nS, nI_BZ, nNB});
947 insertElement({nS, nI_CZ, nI_BZ});
948 insertElement({nS, nI_AZ, nI_CZ});
949 insertElement({nS, nNA, nI_AZ});
950 }
951
952 return true;
953 }
954
955 // Generates geometry for a sharp edge (2 active corners) in 3D (canonical
956 // orientation)
957 template <int Dim = D>
958 std::enable_if_t<Dim == 3, bool> generateCanonicalSharpEdge3D(
959 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
960 int transform, int axis, bool inverted,
961 std::map<hrleIndex, unsigned> *nodes,
962 std::map<hrleIndex, unsigned> *faceNodes,
963 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
964 assert(normalVectorData &&
965 "Normal vector data must be computed for sharp edge generation.");
966
967 // Helper to map global vector to canonical frame
968 auto toCanonical = [&](Vec3D<T> v) {
969 Vec3D<T> res;
970 // 1. Apply reflection (transform)
971 if (transform & 1)
972 v[0] = -v[0];
973 if (transform & 2)
974 v[1] = -v[1];
975 if (transform & 4)
976 v[2] = -v[2];
977 // 2. Apply rotation (axis permutation)
978 // axis=0: XYZ -> XYZ (0,1,2)
979 // axis=1: XYZ -> YZX (1,2,0) - puts Y in X place
980 // axis=2: XYZ -> ZXY (2,0,1) - puts Z in X place
981 if (axis == 0)
982 res = v;
983 else if (axis == 1)
984 res = Vec3D<T>{v[1], v[2], v[0]};
985 else
986 res = Vec3D<T>{v[2], v[0], v[1]};
987 return res;
988 };
989
990 // Helper to map canonical point back to global frame
991 auto fromCanonical = [&](Vec3D<T> v) {
992 Vec3D<T> res;
993 // 1. Inverse rotation
994 if (axis == 0)
995 res = v;
996 else if (axis == 1)
997 res = Vec3D<T>{v[2], v[0], v[1]};
998 else
999 res = Vec3D<T>{v[1], v[2], v[0]};
1000 // 2. Inverse reflection (same as forward for 0/1 swap)
1001 if (transform & 1)
1002 res[0] = T(1.0) - res[0];
1003 if (transform & 2)
1004 res[1] = T(1.0) - res[1];
1005 if (transform & 4)
1006 res[2] = T(1.0) - res[2];
1007 return res;
1008 };
1009
1010 auto getNormal = [&](int idx) {
1011 auto corner = cellIt.getCorner(idx);
1012 if (corner.isDefined()) {
1013 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1014 if (inverted) {
1015 for (int i = 0; i < 3; ++i)
1016 n[i] = -n[i];
1017 }
1018 return n;
1019 }
1020 return Vec3D<T>{};
1021 };
1022
1023 // Calculate face nodes P0 (at x=0) and P1 (at x=1) in canonical frame
1024 Vec3D<T> P[2];
1025 unsigned faceNodeIds[2];
1026 bool nodeExists[2] = {false, false};
1027
1028 // Inverse map indices helper
1029 auto mapIdx = [&](int c) {
1030 // c is canonical index.
1031 // 1. Inverse rotation:
1032 int x = (c & 1), y = (c >> 1) & 1, z = (c >> 2) & 1;
1033 int gx, gy, gz;
1034 if (axis == 0) {
1035 gx = x;
1036 gy = y;
1037 gz = z;
1038 } else if (axis == 1) {
1039 gx = z;
1040 gy = x;
1041 gz = y;
1042 } else {
1043 gx = y;
1044 gy = z;
1045 gz = x;
1046 }
1047
1048 // 2. Inverse reflection
1049 if (transform & 1)
1050 gx = 1 - gx;
1051 if (transform & 2)
1052 gy = 1 - gy;
1053 if (transform & 4)
1054 gz = 1 - gz;
1055
1056 return gx | (gy << 1) | (gz << 2);
1057 };
1058
1059 // First pass: Calculate and check validity without modifying mesh
1060 for (int k = 0; k < 2; ++k) {
1061 // Check if face node already exists
1062 // k=0 -> x=0 face (left). In global, this corresponds to the face
1063 // perpendicular to 'axis' at the 'transform' side. If transform has bit
1064 // 'axis' set, then x=0 in canonical is x=1 in global (relative to cell
1065 // origin). We need to be careful with the face index key. Let's compute
1066 // the global index of the face.
1067 hrleIndex faceIdx = cellIt.getIndices();
1068 // If k=0 (canonical x=0), and transform has bit 'axis' set (flipped),
1069 // this is global x=1 face. If k=1 (canonical x=1), and transform has bit
1070 // 'axis' set, this is global x=0 face.
1071 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1072 if (isHighFace)
1073 faceIdx[axis]++;
1074
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));
1082 }
1083 P[k] = toCanonical(relPos);
1084 } else {
1085 // Calculate P[k]
1086 // For k=0 (x=0): Active corner is 0 (000). Neighbors 2 (010) and 4
1087 // (001). For k=1 (x=1): Active corner is 1 (100). Neighbors 3 (110) and
1088 // 5 (101). We need to map these canonical indices back to global to get
1089 // normals/values.
1090
1091 int c0 = k; // 0 or 1
1092 int c_y = c0 | 2; // Neighbor in Y (canonical)
1093 int c_z = c0 | 4; // Neighbor in Z (canonical)
1094
1095 c0 = mapIdx(c0);
1096 c_y = mapIdx(c_y);
1097 c_z = mapIdx(c_z);
1098
1099 Vec3D<T> n_y = toCanonical(getNormal(c_y));
1100 Vec3D<T> n_z = toCanonical(getNormal(c_z));
1101
1102 if (std::abs(DotProduct(n_y, n_z)) >= 0.8) {
1103 return false;
1104 }
1105
1106 // Solve 2D problem in YZ plane
1107 // Line 1: passes through intersection on edge c0-c_y. Normal n_y
1108 // (projected). Line 2: passes through intersection on edge c0-c_z.
1109 // Normal n_z (projected).
1110 T t_y = getInterp(c0, c_y, cellIt, inverted); // Fraction along Y axis
1111 T t_z = getInterp(c0, c_z, cellIt, inverted); // Fraction along Z axis
1112
1113 // In canonical local frame relative to c0:
1114 // Point on Y-edge: (0, t_y, 0)
1115 // Point on Z-edge: (0, 0, t_z)
1116 // Normals n_y and n_z.
1117 // Solve for (y, z):
1118 // n_y.y * (y - t_y) + n_y.z * (z - 0) = 0
1119 // n_z.y * (y - 0) + n_z.z * (z - t_z) = 0
1120
1121 double det = n_y[1] * n_z[2] - n_y[2] * n_z[1];
1122 if (std::abs(det) < 1e-6)
1123 return false;
1124
1125 double d1 = n_y[1] * t_y;
1126 double d2 = n_z[2] * t_z;
1127
1128 double y = (d1 * n_z[2] - d2 * n_y[2]) / det;
1129 double z = (d2 * n_y[1] - d1 * n_z[1]) / det;
1130
1131 P[k] = Vec3D<T>{(T)k, (T)y, (T)z};
1132
1133 // Check if the point is too far away from the intersection points
1134 // 2 grid deltas = 2.0 in canonical coordinates
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) {
1137 return false;
1138 }
1139 }
1140 }
1141
1142 // Second pass: Commit nodes and stitch
1143 for (int k = 0; k < 2; ++k) {
1144 if (!nodeExists[k]) {
1145 int c0 = k;
1146 hrleIndex faceIdx = cellIt.getIndices();
1147 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1148 if (isHighFace)
1149 faceIdx[axis]++;
1150
1151 Vec3D<T> globalP = fromCanonical(P[k]);
1152 for (int d = 0; d < 3; ++d)
1153 globalP[d] = (globalP[d] + cellIt.getIndices(d));
1154
1155 faceNodeIds[k] = insertNode(globalP);
1156 faceNodes[axis][faceIdx] = faceNodeIds[k];
1157 if (updatePointData && newDataSourceIds)
1158 newDataSourceIds->push_back(
1159 cellIt.getCorner(mapIdx(c0)).getPointId());
1160
1161 stitchToNeighbor(cellIt, axis, isHighFace, faceNodeIds[k], nodes,
1162 valueIt);
1163 }
1164 }
1165
1166 // Calculate edge intersection points in canonical frame
1167 // E02: on edge 0-2 (Y axis). (0, t_02, 0)
1168 // E04: on edge 0-4 (Z axis). (0, 0, t_04)
1169 // E13: on edge 1-3 (Y axis). (1, t_13, 0)
1170 // E15: on edge 1-5 (Z axis). (1, 0, t_15)
1171
1172 // We need to get the node IDs for these. We can use the existing getNode
1173 // helper, but we need to map canonical edges to global edges.
1174 auto getEdgeNode = [&](int c1, int c2) {
1175 // Map canonical corners to global
1176 int g1 = mapIdx(c1);
1177 int g2 = mapIdx(c2);
1178
1179 // Find edge index
1180 int edgeIdx = -1;
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))) {
1186 edgeIdx = e;
1187 break;
1188 }
1189 }
1190 return this->getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1191 };
1192
1193 unsigned n02 = getEdgeNode(0, 2);
1194 unsigned n04 = getEdgeNode(0, 4);
1195 unsigned n13 = getEdgeNode(1, 3);
1196 unsigned n15 = getEdgeNode(1, 5);
1197
1198 // Generate Quads (split into triangles)
1199 // Quad 1 (Y-interface): E02, P0, P1, E13. Normal +Y.
1200 // Winding for +Y normal: E02 -> P0 -> P1 -> E13.
1201 insertElement({n02, faceNodeIds[0], faceNodeIds[1]});
1202 insertElement({n02, faceNodeIds[1], n13});
1203
1204 // Quad 2 (Z-interface): E04, E15, P1, P0. Normal +Z.
1205 // Winding for +Z normal: E04 -> E15 -> P1 -> P0.
1206 insertElement({n04, n15, faceNodeIds[1]});
1207 insertElement({n04, faceNodeIds[1], faceNodeIds[0]});
1208
1209 return true;
1210 }
1211
1212 // Wrapper for 3D sharp edge generation handling different
1213 // rotations/reflections
1214 template <int Dim = D>
1215 std::enable_if_t<Dim == 3, bool> generateSharpEdge3D(
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,
1220 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
1221
1222 bool inverted = false;
1223 int mask = 0;
1224 if (countNeg == 2) {
1225 mask = negMask;
1226 inverted = false;
1227 } else if (countPos == 2) {
1228 mask = posMask;
1229 inverted = true;
1230 } else {
1231 return false;
1232 }
1233
1234 // Check if the two vertices share an edge
1235 int v1 = -1, v2 = -1;
1236 for (int i = 0; i < 8; ++i) {
1237 if ((mask >> i) & 1) {
1238 if (v1 == -1)
1239 v1 = i;
1240 else
1241 v2 = i;
1242 }
1243 }
1244
1245 int diff = v1 ^ v2;
1246 if ((diff & (diff - 1)) != 0)
1247 return false; // Not connected by a single edge
1248
1249 int axis = 0;
1250 while ((diff >> (axis + 1)) > 0)
1251 axis++;
1252
1253 // Transform maps v1 to 0.
1254 int transform = v1;
1255
1256 return generateCanonicalSharpEdge3D(cellIt, transform, axis, inverted,
1257 nodes, faceNodes, newDataSourceIds,
1258 valueIt);
1259 }
1260
1261 // Generates geometry for a sharp corner (1 active corner) in 3D (canonical
1262 // orientation)
1264 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
1265 int transform, bool inverted, std::map<hrleIndex, unsigned> *nodes,
1266 std::map<hrleIndex, unsigned> *faceNodes,
1267 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt,
1268 Vec3D<T> &cornerPos) {
1269
1270 const auto *normalVectorData =
1271 currentLevelSet->getPointData().getVectorData(
1273 assert(normalVectorData &&
1274 "Normal vector data must be available for sharp corner generation");
1275
1276 // Helper to map global vector to canonical frame (just reflection)
1277 auto toCanonical = [&](Vec3D<T> v) {
1278 if (transform & 1)
1279 v[0] = -v[0];
1280 if (transform & 2)
1281 v[1] = -v[1];
1282 if (transform & 4)
1283 v[2] = -v[2];
1284 return v;
1285 };
1286
1287 // Helper to map canonical point back to global frame
1288 auto fromCanonical = [&](Vec3D<T> v) {
1289 if (transform & 1)
1290 v[0] = T(1.0) - v[0];
1291 if (transform & 2)
1292 v[1] = T(1.0) - v[1];
1293 if (transform & 4)
1294 v[2] = T(1.0) - v[2];
1295 return v;
1296 };
1297
1298 auto getNormal = [&](int idx) {
1299 auto corner = cellIt.getCorner(idx);
1300 if (corner.isDefined()) {
1301 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1302 if (inverted) {
1303 for (int i = 0; i < 3; ++i)
1304 n[i] = -n[i];
1305 }
1306 return n;
1307 }
1308 return Vec3D<T>{};
1309 };
1310
1311 Vec3D<T> norm1, norm2, norm3;
1312 double d1, d2, d3;
1313 Vec3D<T> P1, P2, P4;
1314
1315 const int g1 = transform ^ 1;
1316 const int g2 = transform ^ 2;
1317 const int g4 = transform ^ 4;
1318
1319 norm1 = toCanonical(getNormal(g1));
1320 norm2 = toCanonical(getNormal(g2));
1321 norm3 = toCanonical(getNormal(g4));
1322
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);
1326
1327 P1 = {t1, 0, 0};
1328 P2 = {0, t2, 0};
1329 P4 = {0, 0, t4};
1330
1331 d1 = DotProduct(norm1, P1);
1332 d2 = DotProduct(norm2, P2);
1333 d3 = DotProduct(norm3, P4);
1334
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) {
1338 return false;
1339 }
1340
1341 // We solve the system of 3 plane equations:
1342 // norm1 . S = d1
1343 // norm2 . S = d2
1344 // norm3 . S = d3
1345 // This is a linear system Ax = b, where A is the matrix of normals, x is
1346 // the corner position S, and b is the vector of dot products.
1347
1348 // Using Cramer's rule to solve the system
1349
1350 const Vec3D<T> c23 = CrossProduct(norm2, norm3);
1351 const Vec3D<T> c31 = CrossProduct(norm3, norm1);
1352 const Vec3D<T> c12 = CrossProduct(norm1, norm2);
1353
1354 const double det = DotProduct(norm1, c23);
1355
1356 if (std::abs(det) < epsilon)
1357 return false; // Planes are parallel or linearly dependent
1358
1359 const T invDet = 1.0 / det;
1360 Vec3D<T> S;
1361 for (int i = 0; i < 3; ++i) {
1362 S[i] = (d1 * c23[i] + d2 * c31[i] + d3 * c12[i]) * invDet;
1363 }
1364
1365 if (Norm2(S - P1) > 9.0 || Norm2(S - P2) > 9.0 || Norm2(S - P4) > 9.0) {
1366 return false;
1367 }
1368
1369 // Transform corner position back from canonical to global coordinates
1370 Vec3D<T> globalS = fromCanonical(S);
1371 for (int d = 0; d < 3; ++d) {
1372 cornerPos[d] = (globalS[d] + cellIt.getIndices(d));
1373 }
1374
1375 return true;
1376 }
1377
1378 // Wrapper for 3D sharp corner generation handling different
1379 // rotations/reflections
1380 template <int Dim = D>
1381 std::enable_if_t<Dim == 3, bool> generateSharpCorner3D(
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,
1387 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
1388
1389 bool inverted = false;
1390 int transform = -1;
1391 // bool is3vs5 = false;
1392
1393 if (countNeg == 1) {
1394 transform = 0;
1395 while (!((negMask >> transform) & 1))
1396 transform++;
1397 inverted = false;
1398 } else if (countPos == 1) {
1399 transform = 0;
1400 while (!((posMask >> transform) & 1))
1401 transform++;
1402 inverted = true;
1403 }
1404
1405 if (transform == -1)
1406 return false;
1407
1408 Vec3D<T> cornerPos;
1409 if (generateCanonicalSharpCorner3D(cellIt, transform, inverted, nodes,
1410 faceNodes, newDataSourceIds, valueIt,
1411 cornerPos)) {
1412 unsigned nCorner;
1413 hrleIndex cornerIdx = cellIt.getIndices();
1414 cornerIdx += viennahrle::BitMaskToIndex<D>(transform);
1415
1416 auto it = cornerNodes.find(cornerIdx);
1417 if (it != cornerNodes.end()) {
1418 nCorner = it->second;
1419 } else {
1420 nCorner = insertNode(cornerPos);
1421 cornerNodes[cornerIdx] = nCorner;
1422
1423 // Store this sharp corner node for potential use by derived classes
1424 matSharpCornerNodes.push_back({nCorner, cornerPos});
1425
1426 if (updatePointData && newDataSourceIds)
1427 newDataSourceIds->push_back(cellIt.getCorner(transform).getPointId());
1428 }
1429
1430 // Get edge nodes
1431 auto getEdgeNode = [&](int neighbor) {
1432 int v1 = transform;
1433 int v2 = neighbor;
1434 int edgeIdx = -1;
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))) {
1440 edgeIdx = e;
1441 break;
1442 }
1443 }
1444 return this->getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1445 };
1446
1447 unsigned n1 = getEdgeNode(transform ^ 1);
1448 unsigned n2 = getEdgeNode(transform ^ 2);
1449 unsigned n4 = getEdgeNode(transform ^ 4);
1450
1451 // Try to generate cube corner topology
1452 if (normalVectorData) {
1453 // Re-calculate canonical normals and plane constants
1454 auto toCanonical = [&](Vec3D<T> v) {
1455 if (transform & 1)
1456 v[0] = -v[0];
1457 if (transform & 2)
1458 v[1] = -v[1];
1459 if (transform & 4)
1460 v[2] = -v[2];
1461 return v;
1462 };
1463
1464 auto fromCanonical = [&](Vec3D<T> v) {
1465 if (transform & 1)
1466 v[0] = T(1.0) - v[0];
1467 if (transform & 2)
1468 v[1] = T(1.0) - v[1];
1469 if (transform & 4)
1470 v[2] = T(1.0) - v[2];
1471 return v;
1472 };
1473
1474 auto getNormal = [&](int idx) {
1475 auto corner = cellIt.getCorner(idx);
1476 if (corner.isDefined()) {
1477 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1478 if (inverted)
1479 for (int i = 0; i < 3; ++i)
1480 n[i] = -n[i];
1481 return n;
1482 }
1483 return Vec3D<T>{};
1484 };
1485
1486 Vec3D<T> norm1, norm2, norm3;
1487 double d1, d2, d3;
1488
1489 // Identify neighbors
1490 int n1_idx = -1, n2_idx = -1, n3_idx = -1;
1491 int neighbors[] = {transform ^ 1, transform ^ 2, transform ^ 4};
1492
1493 for (int n_idx : neighbors) {
1494 T val = cellIt.getCorner(n_idx).getValue();
1495 // If signs are opposite (product < 0), this is a face neighbor
1496 if (val * cellIt.getCorner(transform).getValue() < 0) {
1497 if (n1_idx == -1)
1498 n1_idx = n_idx;
1499 else
1500 n2_idx = n_idx;
1501 } else {
1502 n3_idx = n_idx;
1503 }
1504 }
1505
1506 // Standard 1-vs-7 case
1507 norm1 = toCanonical(getNormal(transform ^ 1));
1508 norm2 = toCanonical(getNormal(transform ^ 2));
1509 norm3 = toCanonical(getNormal(transform ^ 4));
1510
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);
1514
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)
1519 return false;
1520 res1 = (c1 * b2 - c2 * b1) / det;
1521 res2 = (a1 * c2 - a2 * c1) / det;
1522 return true;
1523 };
1524
1525 Vec3D<T> Fx, Fy, Fz;
1526 Fx[0] = 0.0;
1527 Fy[1] = 0.0;
1528 Fz[2] = 0.0;
1529
1530 bool valid = true;
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]);
1537
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;
1541 };
1542
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)
1547 faceIdx[axis]++;
1548
1549 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
1550 return faceNodes[axis][faceIdx];
1551 }
1552 Vec3D<T> globalPos = fromCanonical(pos);
1553 for (int d = 0; d < 3; ++d)
1554 globalPos[d] = (globalPos[d] + cellIt.getIndices(d));
1555 unsigned nodeId = insertNode(globalPos);
1556 faceNodes[axis][faceIdx] = nodeId;
1557 if (updatePointData && newDataSourceIds)
1558 newDataSourceIds->push_back(
1559 cellIt.getCorner(transform).getPointId());
1560 stitchToNeighbor(cellIt, axis, (transform >> axis) & 1, nodeId,
1561 nodes, valueIt);
1562 return nodeId;
1563 };
1564
1565 unsigned nFx = getFaceNode(0, Fx);
1566 unsigned nFy = getFaceNode(1, Fy);
1567 unsigned nFz = getFaceNode(2, Fz);
1568
1569 // Calculate parity of the transform (number of reflections)
1570 int parity =
1571 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1572
1573 // Flip winding if parity is odd XOR inverted is true.
1574 bool flip = (parity % 2 != 0) ^ inverted;
1575
1576 if (!flip) {
1577 insertElement({nCorner, nFy, n1});
1578 insertElement({nCorner, n1, nFz});
1579 insertElement({nCorner, nFz, n2});
1580 insertElement({nCorner, n2, nFx});
1581 insertElement({nCorner, nFx, n4});
1582 insertElement({nCorner, n4, nFy});
1583 } else {
1584 insertElement({nCorner, n1, nFy});
1585 insertElement({nCorner, nFz, n1});
1586 insertElement({nCorner, n2, nFz});
1587 insertElement({nCorner, nFx, n2});
1588 insertElement({nCorner, n4, nFx});
1589 insertElement({nCorner, nFy, n4});
1590 }
1591 return true;
1592 }
1593 }
1594
1595 // Triangles
1596 // Cycle around normal (1,1,1) in canonical is 1->2->4->1 (x->y->z->x)
1597 // Triangles: (n1, S, n4), (n4, S, n2), (n2, S, n1)
1598 // If inverted, flip winding.
1599
1600 // Calculate parity of the transform (number of reflections)
1601 int parity =
1602 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1603
1604 // Flip winding if parity is odd XOR inverted is true.
1605 bool flip = (parity % 2 != 0) ^ inverted;
1606
1607 if (!flip) {
1608 insertElement({n1, nCorner, n4});
1609 insertElement({n4, nCorner, n2});
1610 insertElement({n2, nCorner, n1});
1611 } else {
1612 insertElement({n1, n4, nCorner});
1613 insertElement({n4, n2, nCorner});
1614 insertElement({n2, n1, nCorner});
1615 }
1616
1617 return true;
1618 }
1619
1620 return false;
1621 }
1622
1623 static inline bool
1624 triangleMisformed(const std::array<unsigned, D> &nodeNumbers) noexcept {
1625 if constexpr (D == 3) {
1626 return nodeNumbers[0] == nodeNumbers[1] ||
1627 nodeNumbers[0] == nodeNumbers[2] ||
1628 nodeNumbers[1] == nodeNumbers[2];
1629 } else {
1630 return nodeNumbers[0] == nodeNumbers[1];
1631 }
1632 }
1633
1634 static inline Vec3D<T> calculateNormal(const Vec3D<T> &nodeA,
1635 const Vec3D<T> &nodeB,
1636 const Vec3D<T> &nodeC) noexcept {
1637 return CrossProduct(nodeB - nodeA, nodeC - nodeA);
1638 }
1639
1640 T getInterp(int p_a, int p_b,
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;
1646 };
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)
1650 return T(0.5);
1651 return v_a / (v_a - v_b);
1652 }
1653
1654 void insertElement(const std::array<unsigned, D> &nodeNumbers) {
1655 if (triangleMisformed(nodeNumbers))
1656 return;
1657
1658 auto elementI3 = [&](const std::array<unsigned, D> &element) -> I3 {
1659 if constexpr (D == 2)
1660 return {(int)element[0], (int)element[1], 0};
1661 else
1662 return {(int)element[0], (int)element[1], (int)element[2]};
1663 };
1664
1665 if (uniqueElements.insert(elementI3(nodeNumbers)).second) {
1666 Vec3D<T> normal;
1667 if constexpr (D == 2) {
1668 normal = Vec3D<T>{
1669 -(mesh->nodes[nodeNumbers[1]][1] - mesh->nodes[nodeNumbers[0]][1]),
1670 mesh->nodes[nodeNumbers[1]][0] - mesh->nodes[nodeNumbers[0]][0],
1671 T(0)};
1672 } else {
1673 normal = calculateNormal(mesh->nodes[nodeNumbers[0]],
1674 mesh->nodes[nodeNumbers[1]],
1675 mesh->nodes[nodeNumbers[2]]);
1676 }
1677
1678 double n2 =
1679 normal[0] * normal[0] + normal[1] * normal[1] + normal[2] * normal[2];
1680 if (n2 > epsilon) {
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++) {
1684 normal[d] *= invn;
1685 }
1686 currentNormals.push_back(normal);
1687 currentMaterials.push_back(static_cast<T>(currentMaterialId));
1688 }
1689 }
1690 }
1691
1692 unsigned insertNode(Vec3D<T> const &pos) {
1693 auto quantize = [&](const Vec3D<T> &p) -> I3 {
1694 const T inv = T(1) / minNodeDistanceFactor;
1695 return {(int)std::llround(p[0] * inv), (int)std::llround(p[1] * inv),
1696 (int)std::llround(p[2] * inv)};
1697 };
1698
1699 int nodeIdx = -1;
1700 if (minNodeDistanceFactor > 0) {
1701 auto q = quantize(pos);
1702 auto it = nodeIdByBin.find(q);
1703 if (it != nodeIdByBin.end())
1704 nodeIdx = it->second;
1705 }
1706
1707 if (nodeIdx >= 0) {
1708 for (int i = 0; i < 3; ++i) {
1709 mesh->nodes[nodeIdx][i] = (mesh->nodes[nodeIdx][i] + pos[i]) * T(0.5);
1710 }
1711 return nodeIdx;
1712 } else {
1713 unsigned newNodeId = mesh->insertNextNode(pos);
1714 if (minNodeDistanceFactor > 0) {
1715 nodeIdByBin.emplace(quantize(pos), newNodeId);
1716 }
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];
1722 }
1723 return newNodeId;
1724 }
1725 }
1726};
1727
1728// add all template specialisations for this class
1729PRECOMPILE_PRECISION_DIMENSION(ToSurfaceMesh)
1730
1731} // namespace viennals
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