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
421 void
422 stitchToNeighbor(viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
423 int axis, bool isHighFace, unsigned faceNodeId,
424 std::map<hrleIndex, unsigned> *nodes,
425 ConstSparseIterator &valueIt) {
426 // Backward stitching: Check if neighbor on this face is "past" and needs
427 // stitching
428 hrleIndex neighborIdx = cellIt.getIndices();
429 if (isHighFace)
430 neighborIdx[axis]++;
431 else
432 neighborIdx[axis]--;
433
434 if (neighborIdx < cellIt.getIndices()) {
435 unsigned nSigns = 0;
436 auto &grid = currentLevelSet->getGrid();
437 for (int i = 0; i < 8; ++i) {
438 hrleIndex cIdx = neighborIdx + viennahrle::BitMaskToIndex<D>(i);
439 if (!grid.isOutsideOfDomain(cIdx)) {
440 valueIt.goToIndices(cIdx);
441 if (valueIt.getValue() >= 0)
442 nSigns |= (1 << i);
443 }
444 }
445
446 if (nSigns != 0 && nSigns != 255) {
447 const int *nTriangles = lsInternal::MarchingCubes::polygonize3d(nSigns);
448 int nFaceD = isHighFace ? 0 : 1;
449
450 for (; nTriangles[0] != -1; nTriangles += 3) {
451 std::vector<unsigned> face_edge_nodes;
452 for (int i = 0; i < 3; ++i) {
453 int edge = nTriangles[i];
454 int c0 = corner0[edge];
455 int c1 = corner1[edge];
456 bool onFace = (((c0 >> axis) & 1) == nFaceD) &&
457 (((c1 >> axis) & 1) == nFaceD);
458 if (onFace) {
459 unsigned p0 = corner0[edge];
460 auto dir = direction[edge];
461 hrleIndex d = neighborIdx + viennahrle::BitMaskToIndex<D>(p0);
462 auto itN = nodes[dir].find(d);
463 if (itN != nodes[dir].end()) {
464 face_edge_nodes.push_back(itN->second);
465 }
466 }
467 }
468 if (face_edge_nodes.size() == 2) {
469 insertElement(std::array<unsigned, 3>{
470 face_edge_nodes[0], face_edge_nodes[1], faceNodeId});
471 }
472 }
473 }
474 }
475 }
476
477 // Solves for a sharp corner position in 2D by intersecting normals
479 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
480 int transform, Vec3D<T> &cornerPos) const {
481 assert(normalVectorData &&
482 "Normal vector data must be computed for sharp corner generation.");
483
484 auto getTransformedGradient = [&](int cornerID) {
485 auto corner = cellIt.getCorner(cornerID ^ transform);
486 if (corner.isDefined()) {
487 auto normal = (*normalVectorData)[corner.getPointId()];
488 normal[2] = 0;
489 if ((transform & 1) != 0)
490 normal[0] = -normal[0];
491 if ((transform & 2) != 0)
492 normal[1] = -normal[1];
493 return normal;
494 }
495 return Vec3D<T>{};
496 };
497
498 Vec3D<T> norm1 = getTransformedGradient(1); // neighbor in x
499 Vec3D<T> norm2 = getTransformedGradient(2); // neighbor in y
500
501 if (std::abs(DotProduct(norm1, norm2)) >= 0.8) {
502 return false;
503 }
504
505 auto calculateNodePos = [&](int edge, Vec3D<T> &pos) {
506 unsigned p0 = corner0[edge];
507 unsigned p1 = corner1[edge];
508 auto dir = direction[edge];
509
510 for (int z = 0; z < D; z++) {
511 if (z != dir) {
512 pos[z] = static_cast<T>(viennahrle::BitMaskToIndex<D>(p0)[z]);
513 } else {
514 const T d0 = cellIt.getCorner(p0 ^ transform).getValue();
515 const T d1 = cellIt.getCorner(p1 ^ transform).getValue();
516 if (d0 == -d1) {
517 pos[z] = T(0.5);
518 } else {
519 if (std::abs(d0) <= std::abs(d1)) {
520 pos[z] = (d0 / (d0 - d1));
521 } else {
522 pos[z] = T(1.0) - (d1 / (d1 - d0));
523 }
524 }
525 pos[z] = std::max(pos[z], epsilon);
526 pos[z] = std::min(pos[z], T(1.0) - epsilon);
527 }
528 }
529 };
530
531 Vec3D<T> pX{}, pY{};
532 calculateNodePos(0, pX); // Edge along x-axis from corner 0
533 calculateNodePos(3, pY); // Edge along y-axis from corner 0
534
535 double d1 = DotProduct(norm1, pX);
536 double d2 = DotProduct(norm2, pY);
537 double det = norm1[0] * norm2[1] - norm1[1] * norm2[0];
538
539 if (std::abs(det) > 1e-6 && std::isfinite(det)) {
540 cornerPos[0] = (d1 * norm2[1] - d2 * norm1[1]) / det;
541 cornerPos[1] = (d2 * norm1[0] - d1 * norm2[0]) / det;
542 } else {
543 for (int i = 0; i < D; ++i)
544 cornerPos[i] = (pX[i] + pY[i]) * T(0.5);
545 }
546
547 return true;
548 }
549
550 // Wrapper for 2D sharp corner generation handling different
551 // rotations/reflections
553 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt, int countNeg,
554 int countPos, int negMask, int posMask,
555 std::map<hrleIndex, unsigned> *nodes,
556 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
557
558 int cornerIdx = -1;
559 if (countNeg == 1) {
560 for (int i = 0; i < 4; ++i)
561 if ((negMask >> i) & 1) {
562 cornerIdx = i;
563 break;
564 }
565 } else if (countPos == 1) {
566 for (int i = 0; i < 4; ++i)
567 if ((posMask >> i) & 1) {
568 cornerIdx = i;
569 break;
570 }
571 }
572
573 if (cornerIdx != -1) {
574 // Check if this corner is also a corner in any neighboring cell
575 bool isSharedCorner = false;
576 auto pIdx =
577 cellIt.getIndices() + viennahrle::BitMaskToIndex<D>(cornerIdx);
578 T pVal = cellIt.getCorner(cornerIdx).getValue();
579 auto &grid = currentLevelSet->getGrid();
580
581 for (int i = 0; i < (1 << D); ++i) {
582 if (i == cornerIdx)
583 continue;
584
585 hrleIndex neighborIndices;
586 for (int k = 0; k < D; ++k)
587 neighborIndices[k] = pIdx[k] - ((i >> k) & 1);
588
589 bool neighborIsCorner = true;
590
591 // Check edge-connected neighbors in the neighbor cell
592 for (int k = 0; k < D; ++k) {
593 int neighborLocal = i ^ (1 << k);
594 auto checkIdx =
595 neighborIndices + viennahrle::BitMaskToIndex<D>(neighborLocal);
596 if (grid.isOutsideOfDomain(checkIdx)) {
597 checkIdx = grid.globalIndices2LocalIndices(checkIdx);
598 }
599 valueIt.goToIndices(checkIdx);
600 T nVal = valueIt.getValue();
601 bool pSign = (pVal >= 0);
602 bool nSign = (nVal >= 0);
603 if (pSign == nSign) {
604 neighborIsCorner = false;
605 break;
606 }
607 }
608
609 if (neighborIsCorner) {
610 isSharedCorner = true;
611 break;
612 }
613 }
614
615 if (!isSharedCorner) {
616 Vec3D<T> cornerPos{};
617 if (generateCanonicalSharpCorner2D(cellIt, cornerIdx, cornerPos)) {
618
619 // inverse transform cornerPos from canonical local to this cell's
620 // local
621 if ((cornerIdx & 1) != 0)
622 cornerPos[0] = T(1.0) - cornerPos[0];
623 if ((cornerIdx & 2) != 0)
624 cornerPos[1] = T(1.0) - cornerPos[1];
625
626 // convert to global coordinates
627 for (int i = 0; i < D; ++i) {
628 cornerPos[i] = (cornerPos[i] + cellIt.getIndices(i));
629 }
630
631 // Determine edges for this corner
632 int edgeX = -1, edgeY = -1;
633 if (cornerIdx == 0) {
634 edgeX = 0;
635 edgeY = 3;
636 } else if (cornerIdx == 1) {
637 edgeX = 0;
638 edgeY = 1;
639 } else if (cornerIdx == 2) {
640 edgeX = 2;
641 edgeY = 3;
642 } else if (cornerIdx == 3) {
643 edgeX = 2;
644 edgeY = 1;
645 }
646
647 unsigned nX = getNode(cellIt, edgeX, nodes, newDataSourceIds);
648 unsigned nY = getNode(cellIt, edgeY, nodes, newDataSourceIds);
649 unsigned nCorner = insertNode(cornerPos);
650
651 // Store this sharp corner node for potential use by derived classes
652 matSharpCornerNodes.push_back({nCorner, cornerPos});
653
654 if (updatePointData && newDataSourceIds) {
655 assert(!newDataSourceIds->empty());
656 newDataSourceIds->push_back(
657 newDataSourceIds->back()); // TODO: improve point data source
658 }
659
660 insertElement({nX, nCorner});
661 insertElement({nCorner, nY});
662 return true;
663 }
664 }
665 }
666 return false;
667 }
668
669 // Generates geometry for an "L-shaped" configuration (3 active corners) in 3D
670 bool
671 generateSharpL3D(viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
672 int countNeg, int countPos, int negMask, int posMask,
673 std::map<hrleIndex, unsigned> *nodes,
674 std::map<hrleIndex, unsigned> *faceNodes,
675 std::vector<unsigned> *newDataSourceIds,
676 ConstSparseIterator &valueIt) {
677
678 bool inverted = false;
679 int mask = 0;
680 if (countNeg == 3) {
681 mask = negMask;
682 inverted = false;
683 } else if (countPos == 3) {
684 mask = posMask;
685 inverted = true;
686 } else {
687 return false;
688 }
689
690 // Check for L-shape: 3 points on a face
691 int C = -1, A = -1, B = -1;
692 for (int i = 0; i < 8; ++i) {
693 if ((mask >> i) & 1) {
694 int neighbors = 0;
695 for (int k = 0; k < 3; ++k) {
696 if ((mask >> (i ^ (1 << k))) & 1)
697 neighbors++;
698 }
699 if (neighbors == 2) {
700 C = i;
701 break;
702 }
703 }
704 }
705
706 if (C == -1)
707 return false;
708
709 // Identify A and B
710 for (int k = 0; k < 3; ++k) {
711 int n = C ^ (1 << k);
712 if ((mask >> n) & 1) {
713 if (A == -1)
714 A = n;
715 else
716 B = n;
717 }
718 }
719
720 assert(normalVectorData &&
721 "Normal vector data must be computed for sharp feature generation.");
722
723 auto getNormal = [&](int idx) {
724 auto corner = cellIt.getCorner(idx);
725 if (corner.isDefined()) {
726 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
727 if (inverted)
728 for (int i = 0; i < 3; ++i)
729 n[i] = -n[i];
730 return n;
731 }
732 return Vec3D<T>{};
733 };
734
735 // Determine axes
736 int axisA = 0;
737 while ((C ^ A) != (1 << axisA))
738 axisA++;
739 int axisB = 0;
740 while ((C ^ B) != (1 << axisB))
741 axisB++;
742 int axisZ = 3 - axisA - axisB;
743
744 // Solve 2D face problem to find sharp point on a face
745 auto calculateFace = [&](int axis, int corner, int n1,
746 int n2) -> std::optional<Vec3D<T>> {
747 hrleIndex faceIdx = cellIt.getIndices();
748 if ((corner >> axis) & 1)
749 faceIdx[axis]++;
750
751 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
752 unsigned nodeId = faceNodes[axis][faceIdx];
753 Vec3D<T> pos = mesh->nodes[nodeId];
754 for (int d = 0; d < 3; ++d)
755 pos[d] = pos[d] - cellIt.getIndices(d);
756 return pos;
757 }
758
759 Vec3D<T> N1 = getNormal(n1);
760 Vec3D<T> N2 = getNormal(n2);
761
762 if (std::abs(DotProduct(N1, N2)) >= 0.8)
763 return std::nullopt;
764
765 int u = (axis + 1) % 3;
766 int v = (axis + 2) % 3;
767
768 int axis1 = 0;
769 while ((corner ^ n1) != (1 << axis1))
770 axis1++;
771 T t1 = getInterp(corner, n1, cellIt, inverted);
772
773 int axis2 = 0;
774 while ((corner ^ n2) != (1 << axis2))
775 axis2++;
776 T t2 = getInterp(corner, n2, cellIt, inverted);
777
778 Vec3D<T> P1;
779 P1[0] = ((corner >> 0) & 1);
780 P1[1] = ((corner >> 1) & 1);
781 P1[2] = ((corner >> 2) & 1);
782 P1[axis1] = ((corner >> axis1) & 1) ? (1.0 - t1) : t1;
783
784 Vec3D<T> P2;
785 P2[0] = ((corner >> 0) & 1);
786 P2[1] = ((corner >> 1) & 1);
787 P2[2] = ((corner >> 2) & 1);
788 P2[axis2] = ((corner >> axis2) & 1) ? (1.0 - t2) : t2;
789
790 double det = N1[u] * N2[v] - N1[v] * N2[u];
791 if (std::abs(det) < 1e-6)
792 return std::nullopt;
793
794 double c1 = N1[u] * P1[u] + N1[v] * P1[v];
795 double c2 = N2[u] * P2[u] + N2[v] * P2[v];
796
797 Vec3D<T> P;
798 P[axis] = P1[axis];
799 P[u] = (c1 * N2[v] - c2 * N1[v]) / det;
800 P[v] = (c2 * N1[u] - c1 * N2[u]) / det;
801
802 if (Norm2(P - P1) > 9.0 || Norm2(P - P2) > 9.0)
803 return std::nullopt;
804
805 return P;
806 };
807
808 auto commitFace = [&](int axis, int corner, Vec3D<T> P) -> unsigned {
809 hrleIndex faceIdx = cellIt.getIndices();
810 if ((corner >> axis) & 1)
811 faceIdx[axis]++;
812
813 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
814 return faceNodes[axis][faceIdx];
815 }
816
817 Vec3D<T> globalP = P;
818 for (int d = 0; d < 3; ++d)
819 globalP[d] = (globalP[d] + cellIt.getIndices(d));
820 unsigned nodeId = insertNode(globalP);
821 faceNodes[axis][faceIdx] = nodeId;
822 if (updatePointData && newDataSourceIds)
823 newDataSourceIds->push_back(cellIt.getCorner(corner).getPointId());
824
825 stitchToNeighbor(cellIt, axis, (corner >> axis) & 1, nodeId, nodes,
826 valueIt);
827
828 return nodeId;
829 };
830
831 int D_corner = A ^ (1 << axisB);
832 int A_z = A ^ (1 << axisZ);
833 int B_z = B ^ (1 << axisZ);
834 int C_z = C ^ (1 << axisZ);
835
836 auto P_A_opt = calculateFace(axisA, A, D_corner, A_z);
837 auto P_B_opt = calculateFace(axisB, B, D_corner, B_z);
838
839 if (!P_A_opt.has_value() || !P_B_opt.has_value())
840 return false;
841
842 auto const &P_A = P_A_opt.value();
843 auto const &P_B = P_B_opt.value();
844
845 // Construct S
846 Vec3D<T> S;
847 S[axisA] = P_B[axisA];
848 S[axisB] = P_A[axisB];
849 S[axisZ] = (P_A[axisZ] + P_B[axisZ]) * 0.5;
850
851 unsigned nS = insertNode([&]() {
852 Vec3D<T> gS = S;
853 for (int i = 0; i < 3; ++i)
854 gS[i] = (gS[i] + cellIt.getIndices(i));
855 return gS;
856 }());
857 if (updatePointData && newDataSourceIds)
858 newDataSourceIds->push_back(cellIt.getCorner(C).getPointId());
859
860 // Calculate S_Face
861 Vec3D<T> S_Face = S;
862 S_Face[axisZ] = static_cast<T>((C >> axisZ) & 1);
863
864 unsigned nS_Face;
865 auto faceIdxZ = cellIt.getIndices();
866 if ((C >> axisZ) & 1)
867 faceIdxZ[axisZ]++;
868
869 if (faceNodes[axisZ].find(faceIdxZ) != faceNodes[axisZ].end()) {
870 nS_Face = faceNodes[axisZ][faceIdxZ];
871 } else {
872 nS_Face = insertNode([&]() {
873 Vec3D<T> gS = S_Face;
874 for (int i = 0; i < 3; ++i)
875 gS[i] = (gS[i] + cellIt.getIndices(i));
876 return gS;
877 }());
878 faceNodes[axisZ][faceIdxZ] = nS_Face;
879 if (updatePointData)
880 newDataSourceIds[0].push_back(cellIt.getCorner(C).getPointId());
881
882 stitchToNeighbor(cellIt, axisZ, (C >> axisZ) & 1, nS_Face, nodes,
883 valueIt);
884 }
885
886 // Get face nodes
887 unsigned nNA = commitFace(axisA, A, P_A);
888 unsigned nNB = commitFace(axisB, B, P_B);
889
890 // Get boundary intersection nodes
891 auto getEdgeNode = [&](int v1, int v2) {
892 int edgeIdx = -1;
893 for (int e = 0; e < 12; ++e) {
894 if ((corner0[e] == static_cast<unsigned>(v1) &&
895 corner1[e] == static_cast<unsigned>(v2)) ||
896 (corner0[e] == static_cast<unsigned>(v2) &&
897 corner1[e] == static_cast<unsigned>(v1))) {
898 edgeIdx = e;
899 break;
900 }
901 }
902 return this->getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
903 };
904
905 unsigned nI_AD = getEdgeNode(A, D_corner);
906 unsigned nI_AZ = getEdgeNode(A, A_z);
907 unsigned nI_BD = getEdgeNode(B, D_corner);
908 unsigned nI_BZ = getEdgeNode(B, B_z);
909 unsigned nI_CZ = getEdgeNode(C, C_z);
910
911 // Winding order
912 int parity = (C & 1) + ((C >> 1) & 1) + ((C >> 2) & 1);
913 bool flip = (parity % 2 != 0) ^ inverted;
914
915 int vA = 0;
916 while ((A ^ C) != (1 << vA))
917 vA++;
918 int vB = 0;
919 while ((B ^ C) != (1 << vB))
920 vB++;
921 bool is_cyclic = ((vA + 1) % 3 == vB);
922 if (!is_cyclic) {
923 std::swap(nNA, nNB);
924 std::swap(nI_AD, nI_BD);
925 std::swap(nI_AZ, nI_BZ);
926 }
927
928 if (!flip) {
929 insertElement({nS, nNA, nI_AD});
930 insertElement({nS, nI_AD, nS_Face});
931 insertElement({nS, nS_Face, nI_BD});
932 insertElement({nS, nI_BD, nNB});
933 insertElement({nS, nNB, nI_BZ});
934 insertElement({nS, nI_BZ, nI_CZ});
935 insertElement({nS, nI_CZ, nI_AZ});
936 insertElement({nS, nI_AZ, nNA});
937 } else {
938 insertElement({nS, nI_AD, nNA});
939 insertElement({nS, nS_Face, nI_AD});
940 insertElement({nS, nI_BD, nS_Face});
941 insertElement({nS, nNB, nI_BD});
942 insertElement({nS, nI_BZ, nNB});
943 insertElement({nS, nI_CZ, nI_BZ});
944 insertElement({nS, nI_AZ, nI_CZ});
945 insertElement({nS, nNA, nI_AZ});
946 }
947
948 return true;
949 }
950
951 // Generates geometry for a sharp edge (2 active corners) in 3D (canonical
952 // orientation)
954 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
955 int transform, int axis, bool inverted,
956 std::map<hrleIndex, unsigned> *nodes,
957 std::map<hrleIndex, unsigned> *faceNodes,
958 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
959 assert(normalVectorData &&
960 "Normal vector data must be computed for sharp edge generation.");
961
962 // Helper to map global vector to canonical frame
963 auto toCanonical = [&](Vec3D<T> v) {
964 Vec3D<T> res;
965 // 1. Apply reflection (transform)
966 if (transform & 1)
967 v[0] = -v[0];
968 if (transform & 2)
969 v[1] = -v[1];
970 if (transform & 4)
971 v[2] = -v[2];
972 // 2. Apply rotation (axis permutation)
973 // axis=0: XYZ -> XYZ (0,1,2)
974 // axis=1: XYZ -> YZX (1,2,0) - puts Y in X place
975 // axis=2: XYZ -> ZXY (2,0,1) - puts Z in X place
976 if (axis == 0)
977 res = v;
978 else if (axis == 1)
979 res = Vec3D<T>{v[1], v[2], v[0]};
980 else
981 res = Vec3D<T>{v[2], v[0], v[1]};
982 return res;
983 };
984
985 // Helper to map canonical point back to global frame
986 auto fromCanonical = [&](Vec3D<T> v) {
987 Vec3D<T> res;
988 // 1. Inverse rotation
989 if (axis == 0)
990 res = v;
991 else if (axis == 1)
992 res = Vec3D<T>{v[2], v[0], v[1]};
993 else
994 res = Vec3D<T>{v[1], v[2], v[0]};
995 // 2. Inverse reflection (same as forward for 0/1 swap)
996 if (transform & 1)
997 res[0] = T(1.0) - res[0];
998 if (transform & 2)
999 res[1] = T(1.0) - res[1];
1000 if (transform & 4)
1001 res[2] = T(1.0) - res[2];
1002 return res;
1003 };
1004
1005 auto getNormal = [&](int idx) {
1006 auto corner = cellIt.getCorner(idx);
1007 if (corner.isDefined()) {
1008 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1009 if (inverted) {
1010 for (int i = 0; i < 3; ++i)
1011 n[i] = -n[i];
1012 }
1013 return n;
1014 }
1015 return Vec3D<T>{};
1016 };
1017
1018 // Calculate face nodes P0 (at x=0) and P1 (at x=1) in canonical frame
1019 Vec3D<T> P[2];
1020 unsigned faceNodeIds[2];
1021 bool nodeExists[2] = {false, false};
1022
1023 // Inverse map indices helper
1024 auto mapIdx = [&](int c) {
1025 // c is canonical index.
1026 // 1. Inverse rotation:
1027 int x = (c & 1), y = (c >> 1) & 1, z = (c >> 2) & 1;
1028 int gx, gy, gz;
1029 if (axis == 0) {
1030 gx = x;
1031 gy = y;
1032 gz = z;
1033 } else if (axis == 1) {
1034 gx = z;
1035 gy = x;
1036 gz = y;
1037 } else {
1038 gx = y;
1039 gy = z;
1040 gz = x;
1041 }
1042
1043 // 2. Inverse reflection
1044 if (transform & 1)
1045 gx = 1 - gx;
1046 if (transform & 2)
1047 gy = 1 - gy;
1048 if (transform & 4)
1049 gz = 1 - gz;
1050
1051 return gx | (gy << 1) | (gz << 2);
1052 };
1053
1054 // First pass: Calculate and check validity without modifying mesh
1055 for (int k = 0; k < 2; ++k) {
1056 // Check if face node already exists
1057 // k=0 -> x=0 face (left). In global, this corresponds to the face
1058 // perpendicular to 'axis' at the 'transform' side. If transform has bit
1059 // 'axis' set, then x=0 in canonical is x=1 in global (relative to cell
1060 // origin). We need to be careful with the face index key. Let's compute
1061 // the global index of the face.
1062 hrleIndex faceIdx = cellIt.getIndices();
1063 // If k=0 (canonical x=0), and transform has bit 'axis' set (flipped),
1064 // this is global x=1 face. If k=1 (canonical x=1), and transform has bit
1065 // 'axis' set, this is global x=0 face.
1066 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1067 if (isHighFace)
1068 faceIdx[axis]++;
1069
1070 auto it = faceNodes[axis].find(faceIdx);
1071 if (it != faceNodes[axis].end()) {
1072 nodeExists[k] = true;
1073 faceNodeIds[k] = it->second;
1074 Vec3D<T> relPos = mesh->nodes[faceNodeIds[k]];
1075 for (int d = 0; d < 3; ++d) {
1076 relPos[d] -= static_cast<T>(cellIt.getIndices(d));
1077 }
1078 P[k] = toCanonical(relPos);
1079 } else {
1080 // Calculate P[k]
1081 // For k=0 (x=0): Active corner is 0 (000). Neighbors 2 (010) and 4
1082 // (001). For k=1 (x=1): Active corner is 1 (100). Neighbors 3 (110) and
1083 // 5 (101). We need to map these canonical indices back to global to get
1084 // normals/values.
1085
1086 int c0 = k; // 0 or 1
1087 int c_y = c0 | 2; // Neighbor in Y (canonical)
1088 int c_z = c0 | 4; // Neighbor in Z (canonical)
1089
1090 c0 = mapIdx(c0);
1091 c_y = mapIdx(c_y);
1092 c_z = mapIdx(c_z);
1093
1094 Vec3D<T> n_y = toCanonical(getNormal(c_y));
1095 Vec3D<T> n_z = toCanonical(getNormal(c_z));
1096
1097 if (std::abs(DotProduct(n_y, n_z)) >= 0.8) {
1098 return false;
1099 }
1100
1101 // Solve 2D problem in YZ plane
1102 // Line 1: passes through intersection on edge c0-c_y. Normal n_y
1103 // (projected). Line 2: passes through intersection on edge c0-c_z.
1104 // Normal n_z (projected).
1105 T t_y = getInterp(c0, c_y, cellIt, inverted); // Fraction along Y axis
1106 T t_z = getInterp(c0, c_z, cellIt, inverted); // Fraction along Z axis
1107
1108 // In canonical local frame relative to c0:
1109 // Point on Y-edge: (0, t_y, 0)
1110 // Point on Z-edge: (0, 0, t_z)
1111 // Normals n_y and n_z.
1112 // Solve for (y, z):
1113 // n_y.y * (y - t_y) + n_y.z * (z - 0) = 0
1114 // n_z.y * (y - 0) + n_z.z * (z - t_z) = 0
1115
1116 double det = n_y[1] * n_z[2] - n_y[2] * n_z[1];
1117 if (std::abs(det) < 1e-6)
1118 return false;
1119
1120 double d1 = n_y[1] * t_y;
1121 double d2 = n_z[2] * t_z;
1122
1123 double y = (d1 * n_z[2] - d2 * n_y[2]) / det;
1124 double z = (d2 * n_y[1] - d1 * n_z[1]) / det;
1125
1126 P[k] = Vec3D<T>{(T)k, (T)y, (T)z};
1127
1128 // Check if the point is too far away from the intersection points
1129 // 2 grid deltas = 2.0 in canonical coordinates
1130 if (Norm2(P[k] - Vec3D<T>{(T)k, t_y, 0}) > 9.0 ||
1131 Norm2(P[k] - Vec3D<T>{(T)k, 0, t_z}) > 9.0) {
1132 return false;
1133 }
1134 }
1135 }
1136
1137 // Second pass: Commit nodes and stitch
1138 for (int k = 0; k < 2; ++k) {
1139 if (!nodeExists[k]) {
1140 int c0 = k;
1141 hrleIndex faceIdx = cellIt.getIndices();
1142 bool isHighFace = (k == 1) ^ ((transform >> axis) & 1);
1143 if (isHighFace)
1144 faceIdx[axis]++;
1145
1146 Vec3D<T> globalP = fromCanonical(P[k]);
1147 for (int d = 0; d < 3; ++d)
1148 globalP[d] = (globalP[d] + cellIt.getIndices(d));
1149
1150 faceNodeIds[k] = insertNode(globalP);
1151 faceNodes[axis][faceIdx] = faceNodeIds[k];
1152 if (updatePointData && newDataSourceIds)
1153 newDataSourceIds->push_back(
1154 cellIt.getCorner(mapIdx(c0)).getPointId());
1155
1156 stitchToNeighbor(cellIt, axis, isHighFace, faceNodeIds[k], nodes,
1157 valueIt);
1158 }
1159 }
1160
1161 // Calculate edge intersection points in canonical frame
1162 // E02: on edge 0-2 (Y axis). (0, t_02, 0)
1163 // E04: on edge 0-4 (Z axis). (0, 0, t_04)
1164 // E13: on edge 1-3 (Y axis). (1, t_13, 0)
1165 // E15: on edge 1-5 (Z axis). (1, 0, t_15)
1166
1167 // We need to get the node IDs for these. We can use the existing getNode
1168 // helper, but we need to map canonical edges to global edges.
1169 auto getEdgeNode = [&](int c1, int c2) {
1170 // Map canonical corners to global
1171 int g1 = mapIdx(c1);
1172 int g2 = mapIdx(c2);
1173
1174 // Find edge index
1175 int edgeIdx = -1;
1176 for (int e = 0; e < 12; ++e) {
1177 if ((corner0[e] == static_cast<unsigned>(g1) &&
1178 corner1[e] == static_cast<unsigned>(g2)) ||
1179 (corner0[e] == static_cast<unsigned>(g2) &&
1180 corner1[e] == static_cast<unsigned>(g1))) {
1181 edgeIdx = e;
1182 break;
1183 }
1184 }
1185 return this->getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1186 };
1187
1188 unsigned n02 = getEdgeNode(0, 2);
1189 unsigned n04 = getEdgeNode(0, 4);
1190 unsigned n13 = getEdgeNode(1, 3);
1191 unsigned n15 = getEdgeNode(1, 5);
1192
1193 // Generate Quads (split into triangles)
1194 // Quad 1 (Y-interface): E02, P0, P1, E13. Normal +Y.
1195 // Winding for +Y normal: E02 -> P0 -> P1 -> E13.
1196 insertElement({n02, faceNodeIds[0], faceNodeIds[1]});
1197 insertElement({n02, faceNodeIds[1], n13});
1198
1199 // Quad 2 (Z-interface): E04, E15, P1, P0. Normal +Z.
1200 // Winding for +Z normal: E04 -> E15 -> P1 -> P0.
1201 insertElement({n04, n15, faceNodeIds[1]});
1202 insertElement({n04, faceNodeIds[1], faceNodeIds[0]});
1203
1204 return true;
1205 }
1206
1207 // Wrapper for 3D sharp edge generation handling different
1208 // rotations/reflections
1210 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt, int countNeg,
1211 int countPos, int negMask, int posMask,
1212 std::map<hrleIndex, unsigned> *nodes,
1213 std::map<hrleIndex, unsigned> *faceNodes,
1214 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
1215
1216 bool inverted = false;
1217 int mask = 0;
1218 if (countNeg == 2) {
1219 mask = negMask;
1220 inverted = false;
1221 } else if (countPos == 2) {
1222 mask = posMask;
1223 inverted = true;
1224 } else {
1225 return false;
1226 }
1227
1228 // Check if the two vertices share an edge
1229 int v1 = -1, v2 = -1;
1230 for (int i = 0; i < 8; ++i) {
1231 if ((mask >> i) & 1) {
1232 if (v1 == -1)
1233 v1 = i;
1234 else
1235 v2 = i;
1236 }
1237 }
1238
1239 int diff = v1 ^ v2;
1240 if ((diff & (diff - 1)) != 0)
1241 return false; // Not connected by a single edge
1242
1243 int axis = 0;
1244 while ((diff >> (axis + 1)) > 0)
1245 axis++;
1246
1247 // Transform maps v1 to 0.
1248 int transform = v1;
1249
1250 return generateCanonicalSharpEdge3D(cellIt, transform, axis, inverted,
1251 nodes, faceNodes, newDataSourceIds,
1252 valueIt);
1253 }
1254
1255 // Generates geometry for a sharp corner (1 active corner) in 3D (canonical
1256 // orientation)
1258 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
1259 int transform, bool inverted, std::map<hrleIndex, unsigned> *nodes,
1260 std::map<hrleIndex, unsigned> *faceNodes,
1261 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt,
1262 Vec3D<T> &cornerPos) {
1263
1264 const auto *normalVectorData =
1265 currentLevelSet->getPointData().getVectorData(
1267 assert(normalVectorData &&
1268 "Normal vector data must be available for sharp corner generation");
1269
1270 // Helper to map global vector to canonical frame (just reflection)
1271 auto toCanonical = [&](Vec3D<T> v) {
1272 if (transform & 1)
1273 v[0] = -v[0];
1274 if (transform & 2)
1275 v[1] = -v[1];
1276 if (transform & 4)
1277 v[2] = -v[2];
1278 return v;
1279 };
1280
1281 // Helper to map canonical point back to global frame
1282 auto fromCanonical = [&](Vec3D<T> v) {
1283 if (transform & 1)
1284 v[0] = T(1.0) - v[0];
1285 if (transform & 2)
1286 v[1] = T(1.0) - v[1];
1287 if (transform & 4)
1288 v[2] = T(1.0) - v[2];
1289 return v;
1290 };
1291
1292 auto getNormal = [&](int idx) {
1293 auto corner = cellIt.getCorner(idx);
1294 if (corner.isDefined()) {
1295 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1296 if (inverted) {
1297 for (int i = 0; i < 3; ++i)
1298 n[i] = -n[i];
1299 }
1300 return n;
1301 }
1302 return Vec3D<T>{};
1303 };
1304
1305 Vec3D<T> norm1, norm2, norm3;
1306 double d1, d2, d3;
1307 Vec3D<T> P1, P2, P4;
1308
1309 const int g1 = transform ^ 1;
1310 const int g2 = transform ^ 2;
1311 const int g4 = transform ^ 4;
1312
1313 norm1 = toCanonical(getNormal(g1));
1314 norm2 = toCanonical(getNormal(g2));
1315 norm3 = toCanonical(getNormal(g4));
1316
1317 const T t1 = getInterp(transform, g1, cellIt, inverted);
1318 const T t2 = getInterp(transform, g2, cellIt, inverted);
1319 const T t4 = getInterp(transform, g4, cellIt, inverted);
1320
1321 P1 = {t1, 0, 0};
1322 P2 = {0, t2, 0};
1323 P4 = {0, 0, t4};
1324
1325 d1 = DotProduct(norm1, P1);
1326 d2 = DotProduct(norm2, P2);
1327 d3 = DotProduct(norm3, P4);
1328
1329 if (std::abs(DotProduct(norm1, norm2)) >= 0.8 ||
1330 std::abs(DotProduct(norm1, norm3)) >= 0.8 ||
1331 std::abs(DotProduct(norm2, norm3)) >= 0.8) {
1332 return false;
1333 }
1334
1335 // We solve the system of 3 plane equations:
1336 // norm1 . S = d1
1337 // norm2 . S = d2
1338 // norm3 . S = d3
1339 // This is a linear system Ax = b, where A is the matrix of normals, x is
1340 // the corner position S, and b is the vector of dot products.
1341
1342 // Using Cramer's rule to solve the system
1343
1344 const Vec3D<T> c23 = CrossProduct(norm2, norm3);
1345 const Vec3D<T> c31 = CrossProduct(norm3, norm1);
1346 const Vec3D<T> c12 = CrossProduct(norm1, norm2);
1347
1348 const double det = DotProduct(norm1, c23);
1349
1350 if (std::abs(det) < epsilon)
1351 return false; // Planes are parallel or linearly dependent
1352
1353 const T invDet = 1.0 / det;
1354 Vec3D<T> S;
1355 for (int i = 0; i < 3; ++i) {
1356 S[i] = (d1 * c23[i] + d2 * c31[i] + d3 * c12[i]) * invDet;
1357 }
1358
1359 if (Norm2(S - P1) > 9.0 || Norm2(S - P2) > 9.0 || Norm2(S - P4) > 9.0) {
1360 return false;
1361 }
1362
1363 // Transform corner position back from canonical to global coordinates
1364 Vec3D<T> globalS = fromCanonical(S);
1365 for (int d = 0; d < 3; ++d) {
1366 cornerPos[d] = (globalS[d] + cellIt.getIndices(d));
1367 }
1368
1369 return true;
1370 }
1371
1372 // Wrapper for 3D sharp corner generation handling different
1373 // rotations/reflections
1375 viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt, int countNeg,
1376 int countPos, int negMask, int posMask,
1377 std::map<hrleIndex, unsigned> *nodes,
1378 std::map<hrleIndex, unsigned> &cornerNodes,
1379 std::map<hrleIndex, unsigned> *faceNodes,
1380 std::vector<unsigned> *newDataSourceIds, ConstSparseIterator &valueIt) {
1381
1382 bool inverted = false;
1383 int transform = -1;
1384 // bool is3vs5 = false;
1385
1386 if (countNeg == 1) {
1387 transform = 0;
1388 while (!((negMask >> transform) & 1))
1389 transform++;
1390 inverted = false;
1391 } else if (countPos == 1) {
1392 transform = 0;
1393 while (!((posMask >> transform) & 1))
1394 transform++;
1395 inverted = true;
1396 }
1397
1398 if (transform == -1)
1399 return false;
1400
1401 Vec3D<T> cornerPos;
1402 if (generateCanonicalSharpCorner3D(cellIt, transform, inverted, nodes,
1403 faceNodes, newDataSourceIds, valueIt,
1404 cornerPos)) {
1405 unsigned nCorner;
1406 hrleIndex cornerIdx = cellIt.getIndices();
1407 cornerIdx += viennahrle::BitMaskToIndex<D>(transform);
1408
1409 auto it = cornerNodes.find(cornerIdx);
1410 if (it != cornerNodes.end()) {
1411 nCorner = it->second;
1412 } else {
1413 nCorner = insertNode(cornerPos);
1414 cornerNodes[cornerIdx] = nCorner;
1415
1416 // Store this sharp corner node for potential use by derived classes
1417 matSharpCornerNodes.push_back({nCorner, cornerPos});
1418
1419 if (updatePointData && newDataSourceIds)
1420 newDataSourceIds->push_back(cellIt.getCorner(transform).getPointId());
1421 }
1422
1423 // Get edge nodes
1424 auto getEdgeNode = [&](int neighbor) {
1425 int v1 = transform;
1426 int v2 = neighbor;
1427 int edgeIdx = -1;
1428 for (int e = 0; e < 12; ++e) {
1429 if ((corner0[e] == static_cast<unsigned>(v1) &&
1430 corner1[e] == static_cast<unsigned>(v2)) ||
1431 (corner0[e] == static_cast<unsigned>(v2) &&
1432 corner1[e] == static_cast<unsigned>(v1))) {
1433 edgeIdx = e;
1434 break;
1435 }
1436 }
1437 return this->getNode(cellIt, edgeIdx, nodes, newDataSourceIds);
1438 };
1439
1440 unsigned n1 = getEdgeNode(transform ^ 1);
1441 unsigned n2 = getEdgeNode(transform ^ 2);
1442 unsigned n4 = getEdgeNode(transform ^ 4);
1443
1444 // Try to generate cube corner topology
1445 if (normalVectorData) {
1446 // Re-calculate canonical normals and plane constants
1447 auto toCanonical = [&](Vec3D<T> v) {
1448 if (transform & 1)
1449 v[0] = -v[0];
1450 if (transform & 2)
1451 v[1] = -v[1];
1452 if (transform & 4)
1453 v[2] = -v[2];
1454 return v;
1455 };
1456
1457 auto fromCanonical = [&](Vec3D<T> v) {
1458 if (transform & 1)
1459 v[0] = T(1.0) - v[0];
1460 if (transform & 2)
1461 v[1] = T(1.0) - v[1];
1462 if (transform & 4)
1463 v[2] = T(1.0) - v[2];
1464 return v;
1465 };
1466
1467 auto getNormal = [&](int idx) {
1468 auto corner = cellIt.getCorner(idx);
1469 if (corner.isDefined()) {
1470 Vec3D<T> n = (*normalVectorData)[corner.getPointId()];
1471 if (inverted)
1472 for (int i = 0; i < 3; ++i)
1473 n[i] = -n[i];
1474 return n;
1475 }
1476 return Vec3D<T>{};
1477 };
1478
1479 Vec3D<T> norm1, norm2, norm3;
1480 double d1, d2, d3;
1481
1482 // Identify neighbors
1483 int n1_idx = -1, n2_idx = -1, n3_idx = -1;
1484 int neighbors[] = {transform ^ 1, transform ^ 2, transform ^ 4};
1485
1486 for (int n_idx : neighbors) {
1487 T val = cellIt.getCorner(n_idx).getValue();
1488 // If signs are opposite (product < 0), this is a face neighbor
1489 if (val * cellIt.getCorner(transform).getValue() < 0) {
1490 if (n1_idx == -1)
1491 n1_idx = n_idx;
1492 else
1493 n2_idx = n_idx;
1494 } else {
1495 n3_idx = n_idx;
1496 }
1497 }
1498
1499 // Standard 1-vs-7 case
1500 norm1 = toCanonical(getNormal(transform ^ 1));
1501 norm2 = toCanonical(getNormal(transform ^ 2));
1502 norm3 = toCanonical(getNormal(transform ^ 4));
1503
1504 d1 = norm1[0] * getInterp(transform, transform ^ 1, cellIt, inverted);
1505 d2 = norm2[1] * getInterp(transform, transform ^ 2, cellIt, inverted);
1506 d3 = norm3[2] * getInterp(transform, transform ^ 4, cellIt, inverted);
1507
1508 auto solve2x2 = [&](double a1, double b1, double c1, double a2,
1509 double b2, double c2, T &res1, T &res2) {
1510 double det = a1 * b2 - a2 * b1;
1511 if (std::abs(det) < 1e-6)
1512 return false;
1513 res1 = (c1 * b2 - c2 * b1) / det;
1514 res2 = (a1 * c2 - a2 * c1) / det;
1515 return true;
1516 };
1517
1518 Vec3D<T> Fx, Fy, Fz;
1519 Fx[0] = 0.0;
1520 Fy[1] = 0.0;
1521 Fz[2] = 0.0;
1522
1523 bool valid = true;
1524 valid &= solve2x2(norm2[1], norm2[2], d2 - norm2[0] * Fx[0], norm3[1],
1525 norm3[2], d3 - norm3[0] * Fx[0], Fx[1], Fx[2]);
1526 valid &= solve2x2(norm1[0], norm1[2], d1 - norm1[1] * Fy[1], norm3[0],
1527 norm3[2], d3 - norm3[1] * Fy[1], Fy[0], Fy[2]);
1528 valid &= solve2x2(norm1[0], norm1[1], d1 - norm1[2] * Fz[2], norm2[0],
1529 norm2[1], d2 - norm2[2] * Fz[2], Fz[0], Fz[1]);
1530
1531 auto checkBounds = [&](Vec3D<T> p) {
1532 return p[0] >= -1e-4 && p[0] <= 1.0 + 1e-4 && p[1] >= -1e-4 &&
1533 p[1] <= 1.0 + 1e-4 && p[2] >= -1e-4 && p[2] <= 1.0 + 1e-4;
1534 };
1535
1536 if (valid && checkBounds(Fx) && checkBounds(Fy) && checkBounds(Fz)) {
1537 auto getFaceNode = [&](int axis, Vec3D<T> pos) {
1538 hrleIndex faceIdx = cellIt.getIndices();
1539 if ((transform >> axis) & 1)
1540 faceIdx[axis]++;
1541
1542 if (faceNodes[axis].find(faceIdx) != faceNodes[axis].end()) {
1543 return faceNodes[axis][faceIdx];
1544 }
1545 Vec3D<T> globalPos = fromCanonical(pos);
1546 for (int d = 0; d < 3; ++d)
1547 globalPos[d] = (globalPos[d] + cellIt.getIndices(d));
1548 unsigned nodeId = insertNode(globalPos);
1549 faceNodes[axis][faceIdx] = nodeId;
1550 if (updatePointData)
1551 newDataSourceIds->push_back(
1552 cellIt.getCorner(transform).getPointId());
1553 stitchToNeighbor(cellIt, axis, (transform >> axis) & 1, nodeId,
1554 nodes, valueIt);
1555 return nodeId;
1556 };
1557
1558 unsigned nFx = getFaceNode(0, Fx);
1559 unsigned nFy = getFaceNode(1, Fy);
1560 unsigned nFz = getFaceNode(2, Fz);
1561
1562 // Calculate parity of the transform (number of reflections)
1563 int parity =
1564 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1565
1566 // Flip winding if parity is odd XOR inverted is true.
1567 bool flip = (parity % 2 != 0) ^ inverted;
1568
1569 if (!flip) {
1570 insertElement({nCorner, nFy, n1});
1571 insertElement({nCorner, n1, nFz});
1572 insertElement({nCorner, nFz, n2});
1573 insertElement({nCorner, n2, nFx});
1574 insertElement({nCorner, nFx, n4});
1575 insertElement({nCorner, n4, nFy});
1576 } else {
1577 insertElement({nCorner, n1, nFy});
1578 insertElement({nCorner, nFz, n1});
1579 insertElement({nCorner, n2, nFz});
1580 insertElement({nCorner, nFx, n2});
1581 insertElement({nCorner, n4, nFx});
1582 insertElement({nCorner, nFy, n4});
1583 }
1584 return true;
1585 }
1586 }
1587
1588 // Triangles
1589 // Cycle around normal (1,1,1) in canonical is 1->2->4->1 (x->y->z->x)
1590 // Triangles: (n1, S, n4), (n4, S, n2), (n2, S, n1)
1591 // If inverted, flip winding.
1592
1593 // Calculate parity of the transform (number of reflections)
1594 int parity =
1595 (transform & 1) + ((transform >> 1) & 1) + ((transform >> 2) & 1);
1596
1597 // Flip winding if parity is odd XOR inverted is true.
1598 bool flip = (parity % 2 != 0) ^ inverted;
1599
1600 if (!flip) {
1601 insertElement({n1, nCorner, n4});
1602 insertElement({n4, nCorner, n2});
1603 insertElement({n2, nCorner, n1});
1604 } else {
1605 insertElement({n1, n4, nCorner});
1606 insertElement({n4, n2, nCorner});
1607 insertElement({n2, n1, nCorner});
1608 }
1609
1610 return true;
1611 }
1612
1613 return false;
1614 }
1615
1616 static inline bool
1617 triangleMisformed(const std::array<unsigned, D> &nodeNumbers) noexcept {
1618 if constexpr (D == 3) {
1619 return nodeNumbers[0] == nodeNumbers[1] ||
1620 nodeNumbers[0] == nodeNumbers[2] ||
1621 nodeNumbers[1] == nodeNumbers[2];
1622 } else {
1623 return nodeNumbers[0] == nodeNumbers[1];
1624 }
1625 }
1626
1627 static inline Vec3D<T> calculateNormal(const Vec3D<T> &nodeA,
1628 const Vec3D<T> &nodeB,
1629 const Vec3D<T> &nodeC) noexcept {
1630 return CrossProduct(nodeB - nodeA, nodeC - nodeA);
1631 }
1632
1633 T getInterp(int p_a, int p_b,
1634 const viennahrle::ConstSparseCellIterator<hrleDomainType> &cellIt,
1635 bool inverted) const {
1636 auto getValue = [&](int idx) {
1637 T val = cellIt.getCorner(idx).getValue();
1638 return inverted ? -val : val;
1639 };
1640 const T v_a = getValue(p_a);
1641 const T v_b = getValue(p_b);
1642 if (std::abs(v_a - v_b) < epsilon)
1643 return T(0.5);
1644 return v_a / (v_a - v_b);
1645 }
1646
1647 void insertElement(const std::array<unsigned, D> &nodeNumbers) {
1648 if (triangleMisformed(nodeNumbers))
1649 return;
1650
1651 auto elementI3 = [&](const std::array<unsigned, D> &element) -> I3 {
1652 if constexpr (D == 2)
1653 return {(int)element[0], (int)element[1], 0};
1654 else
1655 return {(int)element[0], (int)element[1], (int)element[2]};
1656 };
1657
1658 if (uniqueElements.insert(elementI3(nodeNumbers)).second) {
1659 Vec3D<T> normal;
1660 if constexpr (D == 2) {
1661 normal = Vec3D<T>{
1662 -(mesh->nodes[nodeNumbers[1]][1] - mesh->nodes[nodeNumbers[0]][1]),
1663 mesh->nodes[nodeNumbers[1]][0] - mesh->nodes[nodeNumbers[0]][0],
1664 T(0)};
1665 } else {
1666 normal = calculateNormal(mesh->nodes[nodeNumbers[0]],
1667 mesh->nodes[nodeNumbers[1]],
1668 mesh->nodes[nodeNumbers[2]]);
1669 }
1670
1671 double n2 =
1672 normal[0] * normal[0] + normal[1] * normal[1] + normal[2] * normal[2];
1673 if (n2 > epsilon) {
1674 mesh->insertNextElement(nodeNumbers);
1675 T invn = static_cast<T>(1.) / std::sqrt(static_cast<T>(n2));
1676 for (int d = 0; d < D; d++) {
1677 normal[d] *= invn;
1678 }
1679 currentNormals.push_back(normal);
1680 currentMaterials.push_back(static_cast<T>(currentMaterialId));
1681 }
1682 }
1683 }
1684
1685 unsigned insertNode(Vec3D<T> const &pos) {
1686 auto quantize = [&](const Vec3D<T> &p) -> I3 {
1687 const T inv = T(1) / minNodeDistanceFactor;
1688 return {(int)std::llround(p[0] * inv), (int)std::llround(p[1] * inv),
1689 (int)std::llround(p[2] * inv)};
1690 };
1691
1692 int nodeIdx = -1;
1693 if (minNodeDistanceFactor > 0) {
1694 auto q = quantize(pos);
1695 auto it = nodeIdByBin.find(q);
1696 if (it != nodeIdByBin.end())
1697 nodeIdx = it->second;
1698 }
1699
1700 if (nodeIdx >= 0) {
1701 for (int i = 0; i < 3; ++i) {
1702 mesh->nodes[nodeIdx][i] = (mesh->nodes[nodeIdx][i] + pos[i]) * T(0.5);
1703 }
1704 return nodeIdx;
1705 } else {
1706 unsigned newNodeId = mesh->insertNextNode(pos);
1707 if (minNodeDistanceFactor > 0) {
1708 nodeIdByBin.emplace(quantize(pos), newNodeId);
1709 }
1710 for (int a = 0; a < D; a++) {
1711 if (pos[a] < mesh->minimumExtent[a])
1712 mesh->minimumExtent[a] = pos[a];
1713 if (pos[a] > mesh->maximumExtent[a])
1714 mesh->maximumExtent[a] = pos[a];
1715 }
1716 return newNodeId;
1717 }
1718 }
1719};
1720
1721// add all template specialisations for this class
1722PRECOMPILE_PRECISION_DIMENSION(ToSurfaceMesh)
1723
1724} // 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
static constexpr unsigned int direction[12]
Definition lsToSurfaceMesh.hpp:77
bool generateSharpEdge3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:1209
void setMesh(SmartPointer< Mesh< T > > passedMesh)
Definition lsToSurfaceMesh.hpp:95
SmartPointer< lsDomainType > currentLevelSet
Definition lsToSurfaceMesh.hpp:67
bool generateSharpCorner2D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:552
bool generateSharpCorners
Definition lsToSurfaceMesh.hpp:38
bool generateCanonicalSharpEdge3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int transform, int axis, bool inverted, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:953
std::vector< SmartPointer< lsDomainType > > levelSets
Definition lsToSurfaceMesh.hpp:32
typename lsDomainType::DomainType hrleDomainType
Definition lsToSurfaceMesh.hpp:28
static bool triangleMisformed(const std::array< unsigned, D > &nodeNumbers) noexcept
Definition lsToSurfaceMesh.hpp:1617
std::unordered_map< I3, unsigned, I3Hash > nodeIdByBin
Definition lsToSurfaceMesh.hpp:61
const T epsilon
Definition lsToSurfaceMesh.hpp:35
viennals::Domain< T, D > lsDomainType
Definition lsToSurfaceMesh.hpp:27
viennahrle::ConstSparseIterator< hrleDomainType > ConstSparseIterator
Definition lsToSurfaceMesh.hpp:30
Vec3D< T > computeNodePosition(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int edge)
Definition lsToSurfaceMesh.hpp:356
bool generateCanonicalSharpCorner2D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int transform, Vec3D< T > &cornerPos) const
Definition lsToSurfaceMesh.hpp:478
bool generateSharpCorner3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > &cornerNodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:1374
void stitchToNeighbor(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int axis, bool isHighFace, unsigned faceNodeId, std::map< hrleIndex, unsigned > *nodes, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:422
std::unordered_set< I3, I3Hash > uniqueElements
Definition lsToSurfaceMesh.hpp:62
std::vector< std::pair< unsigned, Vec3D< T > > > matSharpCornerNodes
Definition lsToSurfaceMesh.hpp:71
static constexpr unsigned int corner1[12]
Definition lsToSurfaceMesh.hpp:75
static Vec3D< T > calculateNormal(const Vec3D< T > &nodeA, const Vec3D< T > &nodeB, const Vec3D< T > &nodeC) noexcept
Definition lsToSurfaceMesh.hpp:1627
const double minNodeDistanceFactor
Definition lsToSurfaceMesh.hpp:36
void insertElement(const std::array< unsigned, D > &nodeNumbers)
Definition lsToSurfaceMesh.hpp:1647
void setUpdatePointData(bool update)
Definition lsToSurfaceMesh.hpp:97
ToSurfaceMesh(double mnd=0.05, double eps=1e-12)
Definition lsToSurfaceMesh.hpp:81
double currentGridDelta
Definition lsToSurfaceMesh.hpp:65
unsigned insertNode(Vec3D< T > const &pos)
Definition lsToSurfaceMesh.hpp:1685
viennahrle::Index< D > hrleIndex
Definition lsToSurfaceMesh.hpp:29
bool updatePointData
Definition lsToSurfaceMesh.hpp:37
bool generateCanonicalSharpCorner3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int transform, bool inverted, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt, Vec3D< T > &cornerPos)
Definition lsToSurfaceMesh.hpp:1257
T getInterp(int p_a, int p_b, const viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, bool inverted) const
Definition lsToSurfaceMesh.hpp:1633
bool generateSharpL3D(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int countNeg, int countPos, int negMask, int posMask, std::map< hrleIndex, unsigned > *nodes, std::map< hrleIndex, unsigned > *faceNodes, std::vector< unsigned > *newDataSourceIds, ConstSparseIterator &valueIt)
Definition lsToSurfaceMesh.hpp:671
PointData< T >::VectorDataType * normalVectorData
Definition lsToSurfaceMesh.hpp:68
std::vector< Vec3D< T > > currentNormals
Definition lsToSurfaceMesh.hpp:63
void scaleMesh()
Definition lsToSurfaceMesh.hpp:342
std::vector< T > currentMaterials
Definition lsToSurfaceMesh.hpp:64
void setSharpCorners(bool check)
Definition lsToSurfaceMesh.hpp:99
virtual void apply()
Definition lsToSurfaceMesh.hpp:101
T currentMaterialId
Definition lsToSurfaceMesh.hpp:66
void setLevelSet(SmartPointer< lsDomainType > passedlsDomain)
Definition lsToSurfaceMesh.hpp:90
unsigned getNode(viennahrle::ConstSparseCellIterator< hrleDomainType > &cellIt, int edge, std::map< hrleIndex, unsigned > *nodes, std::vector< unsigned > *newDataSourceIds)
Definition lsToSurfaceMesh.hpp:386
static constexpr unsigned int corner0[12]
Definition lsToSurfaceMesh.hpp:73
#define PRECOMPILE_PRECISION_DIMENSION(className)
Definition lsPreCompileMacros.hpp:24
Definition lsAdvect.hpp:41
@ ONE_SIDED_MIN_MOD
Definition lsCalculateNormalVectors.hpp:22
Definition lsToSurfaceMesh.hpp:47
size_t operator()(const I3 &k) const
Definition lsToSurfaceMesh.hpp:48
Definition lsToSurfaceMesh.hpp:40
bool operator==(const I3 &o) const
Definition lsToSurfaceMesh.hpp:42
int x
Definition lsToSurfaceMesh.hpp:41
int y
Definition lsToSurfaceMesh.hpp:41
int z
Definition lsToSurfaceMesh.hpp:41