ViennaLS
Loading...
Searching...
No Matches
lsToMultiSurfaceMesh.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <lsMaterialMap.hpp>
4#include <lsToSurfaceMesh.hpp>
5
6namespace viennals {
7
8using namespace viennacore;
9
10template <class NumericType, int D>
11class ToMultiSurfaceMesh : public ToSurfaceMesh<NumericType, D> {
12 using lsDomainType = viennals::Domain<NumericType, D>;
13 using hrleDomainType = typename lsDomainType::DomainType;
14 using hrleIndex = viennahrle::Index<D>;
15
29
30 SmartPointer<MaterialMap> materialMap = nullptr;
31
32 // Structure to track sharp corner nodes from each material
33 struct SharpCornerNode {
34 unsigned nodeId;
35 Vec3D<NumericType> position;
36 };
37
38 // Find the closest sharp corner node from materials below l, subject to a
39 // grid-level layer thickness check. Returns the node ID or -1 if none found.
40 //
41 // activeGridIdx : global grid index of the "active" (single-minority sign)
42 // corner of the cell in which material l's sharp corner was generated.
43 // lsLAtActive : LS value of material l at that grid corner.
44 //
45 // For conformal deposition of material l (thickness t) on material m:
46 // LS_l = LS_m - t everywhere, so
47 // delta = LS_m(activeGridIdx) - LS_l(activeGridIdx) ≈ t.
48 // Skipping snap when delta > layerThreshold prevents merging corners that
49 // belong to physically distinct surfaces separated by a real layer.
50 int findNearestLowerCorner(
51 unsigned l, const Vec3D<NumericType> &pos,
52 const std::unordered_map<unsigned, std::vector<SharpCornerNode>>
53 &sharpCornerNodes,
54 const hrleIndex &activeGridIdx, NumericType lsLAtActive) {
55 int snapId = -1;
57 // Minimum layer thickness (in grid units) that suppresses snapping.
58 const NumericType layerThreshold = NumericType(this->minNodeDistanceFactor);
59 for (unsigned m = 0; m < l; ++m) {
60 auto it = sharpCornerNodes.find(m);
61 if (it == sharpCornerNodes.end())
62 continue;
63 // Evaluate LS_m at the active grid corner of material l's corner cell.
64 // delta approximates the local layer thickness between material m and l.
65 viennahrle::ConstSparseIterator<hrleDomainType> mIt(
66 levelSets[m]->getDomain());
67 mIt.goToIndices(activeGridIdx);
68 const NumericType lsMAtActive = mIt.getValue();
69 // Only apply the layer-thickness check when the active corner is inside
70 // BOTH materials (LS_m < 0 and LS_l < 0). If LS_m > 0, the active corner
71 // lies outside material m (e.g. inside a trench cavity), meaning the two
72 // surfaces meet at a genuine junction — not a thin-layer overlap — and
73 // the snap should be allowed regardless of delta.
74 if (lsMAtActive < NumericType(0) &&
75 lsMAtActive - lsLAtActive > layerThreshold)
76 continue; // thin conformal layer present — do not snap
77 for (const auto &scn : it->second) {
78 NumericType dist2 = NumericType(0);
79 for (int i = 0; i < D; ++i) {
80 NumericType d = scn.position[i] - pos[i];
81 dist2 += d * d;
82 }
83 if (dist2 < minDist2) {
84 minDist2 = dist2;
85 snapId = scn.nodeId;
86 }
87 }
88 }
89 return snapId;
90 }
91
92 // In 2D, check if cornerPos lies on any edge of a previous material's mesh.
93 // If so, split that edge to include cornerNodeId.
94 void splitPreviousMaterialEdge(unsigned l, unsigned cornerNodeId,
95 const Vec3D<NumericType> &cornerPos,
96 size_t meshSizeBefore) {
97 if constexpr (D != 2)
98 return;
99 if (l == 0)
100 return;
101
102 for (size_t lineIdx = 0; lineIdx < mesh->lines.size(); ++lineIdx) {
103 unsigned mat = currentMaterials[lineIdx];
104 if (mat >= l)
105 continue;
106
107 auto &line = mesh->lines[lineIdx];
108 unsigned node0 = line[0];
109 unsigned node1 = line[1];
110
111 auto const &pos0 = mesh->nodes[node0];
112 auto const &pos1 = mesh->nodes[node1];
113 auto lineVec = pos1 - pos0;
114
115 NumericType lineLen2 = NumericType(0);
116 NumericType dot = NumericType(0);
117 for (int i = 0; i < D; ++i) {
118 lineLen2 += lineVec[i] * lineVec[i];
119 dot += (cornerPos[i] - pos0[i]) * lineVec[i];
120 }
121
122 if (lineLen2 < epsilon)
123 continue;
124
125 NumericType t = dot / lineLen2;
126 if (t <= NumericType(0) || t >= NumericType(1))
127 continue;
128
129 // Distance from corner to closest point on line
130 NumericType dist2 = NumericType(0);
131 for (int i = 0; i < D; ++i) {
132 NumericType d = cornerPos[i] - (pos0[i] + t * lineVec[i]);
133 dist2 += d * d;
134 }
135 if (dist2 >= NumericType(this->minNodeDistanceFactor))
136 continue;
137
138 // Split the edge unconditionally to ensure the bottom material remains
139 // closed
140 uniqueElements.erase({(int)node0, (int)node1, 0});
141 uniqueElements.insert({(int)node0, (int)cornerNodeId, 0});
142 uniqueElements.insert({(int)cornerNodeId, (int)node1, 0});
143
144 line[1] = cornerNodeId;
145 mesh->lines.push_back({cornerNodeId, node1});
146 currentMaterials.push_back(mat);
147 Vec3D<NumericType> normal = currentNormals[lineIdx];
148 currentNormals.push_back(normal);
149
150 // Remove any duplicate segments from the current material (l)
151 for (size_t k = meshSizeBefore; k < mesh->lines.size();) {
152 if (currentMaterials[k] != static_cast<NumericType>(l)) {
153 ++k;
154 continue;
155 }
156
157 const auto &rl = mesh->lines[k];
158 bool isSeg1 = (rl[0] == node0 && rl[1] == cornerNodeId) ||
159 (rl[0] == cornerNodeId && rl[1] == node0);
160 bool isSeg2 = (rl[0] == cornerNodeId && rl[1] == node1) ||
161 (rl[0] == node1 && rl[1] == cornerNodeId);
162
163 if (isSeg1 || isSeg2) {
164 mesh->lines.erase(mesh->lines.begin() + k);
165 currentMaterials.erase(currentMaterials.begin() + k);
166 currentNormals.erase(currentNormals.begin() + k);
167 } else {
168 ++k;
169 }
170 }
171
172 break; // Corner can only be on one edge
173 }
174 }
175
176 bool getSegmentIntersection(const Vec3D<NumericType> &p1,
177 const Vec3D<NumericType> &p2,
178 const Vec3D<NumericType> &p3,
179 const Vec3D<NumericType> &p4,
180 Vec3D<NumericType> &intersection,
181 NumericType &t) const {
182 NumericType x1 = p1[0], y1 = p1[1];
183 NumericType x2 = p2[0], y2 = p2[1];
184 NumericType x3 = p3[0], y3 = p3[1];
185 NumericType x4 = p4[0], y4 = p4[1];
186
187 NumericType denom = (y4 - y3) * (x2 - x1) - (x4 - x3) * (y2 - y1);
188 if (std::abs(denom) < epsilon)
189 return false;
190
191 NumericType ua = ((x4 - x3) * (y1 - y3) - (y4 - y3) * (x1 - x3)) / denom;
192 NumericType ub = ((x2 - x1) * (y1 - y3) - (y2 - y1) * (x1 - x3)) / denom;
193
194 if (ua > epsilon && ua < 1.0 - epsilon && ub > epsilon &&
195 ub < 1.0 - epsilon) {
196 intersection[0] = x1 + ua * (x2 - x1);
197 intersection[1] = y1 + ua * (y2 - y1);
198 if constexpr (D == 3)
199 intersection[2] = 0;
200 t = ua;
201 return true;
202 }
203 return false;
204 }
205
206 // In 2D, check if the polymer edge at pIdx crosses BOTH corner edges of any
207 // lower material corner node that lies in the same cell (cellIndices).
208 // When both crossings are found simultaneously:
209 // - Split each mask corner edge at its intersection node.
210 // - Keep Part A (pNode0→intNode1) and Part C (intNode2→pNode1) of the
211 // polymer edge.
212 // - DISCARD Part B (intNode1→intNode2), which lies fictitiously inside the
213 // lower material.
214 void handleCellCornerCrossings(
215 size_t pIdx,
216 const std::unordered_map<unsigned, std::vector<SharpCornerNode>>
217 &sharpCornerNodes,
218 unsigned l, const hrleIndex &cellIndices) {
219 if constexpr (D != 2)
220 return;
221
222 unsigned pNode0 = mesh->lines[pIdx][0];
223 unsigned pNode1 = mesh->lines[pIdx][1];
224 Vec3D<NumericType> pp0 = mesh->nodes[pNode0];
225 Vec3D<NumericType> pp1 = mesh->nodes[pNode1];
226
227 for (unsigned m = 0; m < l; ++m) {
228 auto cit = sharpCornerNodes.find(m);
229 if (cit == sharpCornerNodes.end())
230 continue;
231
232 for (const auto &mscn : cit->second) {
233 // Cell-local check: the lower material's corner must be in this cell
234 bool inCell = true;
235 for (int i = 0; i < D; ++i) {
236 if (mscn.position[i] < NumericType(cellIndices[i]) - epsilon ||
237 mscn.position[i] > NumericType(cellIndices[i] + 1) + epsilon) {
238 inCell = false;
239 break;
240 }
241 }
242 if (!inCell)
243 continue;
244
245 // Find the two corner edges of material m incident to this corner node
246 size_t edge1Idx = SIZE_MAX, edge2Idx = SIZE_MAX;
247 for (size_t eIdx = 0; eIdx < pIdx; ++eIdx) {
248 if (static_cast<unsigned>(currentMaterials[eIdx]) != m)
249 continue;
250 if (mesh->lines[eIdx][0] != mscn.nodeId &&
251 mesh->lines[eIdx][1] != mscn.nodeId)
252 continue;
253 if (edge1Idx == SIZE_MAX)
254 edge1Idx = eIdx;
255 else if (edge2Idx == SIZE_MAX) {
256 edge2Idx = eIdx;
257 break;
258 }
259 }
260 if (edge1Idx == SIZE_MAX || edge2Idx == SIZE_MAX)
261 continue;
262
263 // Check if the polymer edge crosses BOTH corner edges
264 Vec3D<NumericType> int1{}, int2{};
265 NumericType t1 = 0, t2 = 0;
266 bool cross1 = getSegmentIntersection(
267 pp0, pp1, mesh->nodes[mesh->lines[edge1Idx][0]],
268 mesh->nodes[mesh->lines[edge1Idx][1]], int1, t1);
269 bool cross2 = getSegmentIntersection(
270 pp0, pp1, mesh->nodes[mesh->lines[edge2Idx][0]],
271 mesh->nodes[mesh->lines[edge2Idx][1]], int2, t2);
272
273 if (!cross1 || !cross2)
274 continue;
275
276 // Sort by parameter along polymer edge so t1 < t2
277 if (t1 > t2) {
278 std::swap(t1, t2);
279 std::swap(int1, int2);
280 std::swap(edge1Idx, edge2Idx);
281 }
282
283 // Insert intersection nodes
284 unsigned intNode1 = this->insertNode(int1);
285 unsigned intNode2 = this->insertNode(int2);
286
287 // Save metadata before any array modifications
288 NumericType matM1 = currentMaterials[edge1Idx];
289 Vec3D<NumericType> normalM1 = currentNormals[edge1Idx];
290 NumericType matM2 = currentMaterials[edge2Idx];
291 Vec3D<NumericType> normalM2 = currentNormals[edge2Idx];
292 NumericType matL = currentMaterials[pIdx];
293 Vec3D<NumericType> normalL = currentNormals[pIdx];
294 unsigned e1n0 = mesh->lines[edge1Idx][0];
295 unsigned e1n1 = mesh->lines[edge1Idx][1];
296 unsigned e2n0 = mesh->lines[edge2Idx][0];
297 unsigned e2n1 = mesh->lines[edge2Idx][1];
298
299 // Split mask edge1 at intNode1
300 uniqueElements.erase({(int)e1n0, (int)e1n1, 0});
301 mesh->lines[edge1Idx][1] = intNode1;
302 uniqueElements.insert({(int)e1n0, (int)intNode1, 0});
303 mesh->lines.push_back({intNode1, e1n1});
304 currentMaterials.push_back(matM1);
305 currentNormals.push_back(normalM1);
306 uniqueElements.insert({(int)intNode1, (int)e1n1, 0});
307
308 // Split mask edge2 at intNode2
309 uniqueElements.erase({(int)e2n0, (int)e2n1, 0});
310 mesh->lines[edge2Idx][1] = intNode2;
311 uniqueElements.insert({(int)e2n0, (int)intNode2, 0});
312 mesh->lines.push_back({intNode2, e2n1});
313 currentMaterials.push_back(matM2);
314 currentNormals.push_back(normalM2);
315 uniqueElements.insert({(int)intNode2, (int)e2n1, 0});
316
317 // Polymer edge:
318 // Part A (pNode0→intNode1) — kept
319 // Part B (intNode1→intNode2) — DISCARDED (fictitiously inside mask)
320 // Part C (intNode2→pNode1) — kept
321 uniqueElements.erase({(int)pNode0, (int)pNode1, 0});
322 mesh->lines[pIdx][1] = intNode1;
323 uniqueElements.insert({(int)pNode0, (int)intNode1, 0});
324 mesh->lines.push_back({intNode2, pNode1});
325 currentMaterials.push_back(matL);
326 currentNormals.push_back(normalL);
327 uniqueElements.insert({(int)intNode2, (int)pNode1, 0});
328
329 return; // Only one lower-material corner can be in this cell
330 }
331 }
332 }
333
334 // In 2D, for a just-inserted regular MC edge at lineIdx, iterates
335 // lower-material sharp corner nodes once and performs two checks per node:
336 // 1. Point-on-edge: if a corner lies on the edge, split it to close the
337 // shared boundary between adjacent materials.
338 // 2. Double-crossing (cell-local): if the edge crosses both arms of the
339 // corner's L-shape, discard the middle segment that lies fictitiously
340 // inside the lower material (see handleCellCornerCrossings for the same
341 // logic applied to polymer sharp-corner edges).
342 // Returns after the first hit of either kind.
343 void handleEdgeCornerInteractions(
344 size_t lineIdx,
345 const std::unordered_map<unsigned, std::vector<SharpCornerNode>>
346 &sharpCornerNodes,
347 unsigned l, const hrleIndex &cellIndices) {
348 if constexpr (D != 2)
349 return;
350 if (l == 0 || lineIdx >= mesh->lines.size())
351 return;
352
353 unsigned node0 = mesh->lines[lineIdx][0];
354 unsigned node1 = mesh->lines[lineIdx][1];
355 // Copy positions — insertNode may reallocate mesh->nodes
356 Vec3D<NumericType> pos0 = mesh->nodes[node0];
357 Vec3D<NumericType> pos1 = mesh->nodes[node1];
358 auto lineVec = pos1 - pos0;
359
360 NumericType lineLen2 = NumericType(0);
361 for (int i = 0; i < D; ++i)
362 lineLen2 += lineVec[i] * lineVec[i];
363 if (lineLen2 < epsilon)
364 return;
365
366 for (unsigned m = 0; m < l; ++m) {
367 auto it = sharpCornerNodes.find(m);
368 if (it == sharpCornerNodes.end())
369 continue;
370
371 for (const auto &scn : it->second) {
372
373 // --- Check 1: point-on-edge (shared boundary) ---
374 if (scn.nodeId != node0 && scn.nodeId != node1) {
375 NumericType dot = NumericType(0);
376 for (int i = 0; i < D; ++i)
377 dot += (scn.position[i] - pos0[i]) * lineVec[i];
378 NumericType t = dot / lineLen2;
379 if (t > NumericType(0) && t < NumericType(1)) {
380 NumericType dist2 = NumericType(0);
381 for (int i = 0; i < D; ++i) {
382 NumericType d = scn.position[i] - (pos0[i] + t * lineVec[i]);
383 dist2 += d * d;
384 }
385 if (dist2 < NumericType(this->minNodeDistanceFactor)) {
386 unsigned cId = scn.nodeId;
387 bool keep0 = !uniqueElements.count({(int)cId, (int)node0, 0}) &&
388 !uniqueElements.count({(int)node0, (int)cId, 0});
389 bool keep1 = !uniqueElements.count({(int)node1, (int)cId, 0}) &&
390 !uniqueElements.count({(int)cId, (int)node1, 0});
391 uniqueElements.erase({(int)node0, (int)node1, 0});
392 auto &line = mesh->lines[lineIdx];
393 if (keep0 && keep1) {
394 line[1] = cId;
395 mesh->lines.push_back({cId, node1});
396 currentMaterials.push_back(currentMaterials[lineIdx]);
397 currentNormals.push_back(currentNormals[lineIdx]);
398 uniqueElements.insert({(int)node0, (int)cId, 0});
399 uniqueElements.insert({(int)cId, (int)node1, 0});
400 } else if (keep0) {
401 line[1] = cId;
402 uniqueElements.insert({(int)node0, (int)cId, 0});
403 } else if (keep1) {
404 line[0] = cId;
405 uniqueElements.insert({(int)cId, (int)node1, 0});
406 } else {
407 mesh->lines.erase(mesh->lines.begin() + lineIdx);
408 currentMaterials.erase(currentMaterials.begin() + lineIdx);
409 currentNormals.erase(currentNormals.begin() + lineIdx);
410 }
411 return;
412 }
413 }
414 }
415
416 // --- Check 2: double-crossing (cell-local) ---
417 bool inCell = true;
418 for (int i = 0; i < D; ++i) {
419 if (scn.position[i] < NumericType(cellIndices[i]) - epsilon ||
420 scn.position[i] > NumericType(cellIndices[i] + 1) + epsilon) {
421 inCell = false;
422 break;
423 }
424 }
425 if (!inCell)
426 continue;
427
428 size_t edge1Idx = SIZE_MAX, edge2Idx = SIZE_MAX;
429 for (size_t eIdx = 0; eIdx < lineIdx; ++eIdx) {
430 if (static_cast<unsigned>(currentMaterials[eIdx]) != m)
431 continue;
432 if (mesh->lines[eIdx][0] != scn.nodeId &&
433 mesh->lines[eIdx][1] != scn.nodeId)
434 continue;
435 if (edge1Idx == SIZE_MAX)
436 edge1Idx = eIdx;
437 else if (edge2Idx == SIZE_MAX) {
438 edge2Idx = eIdx;
439 break;
440 }
441 }
442 if (edge1Idx == SIZE_MAX || edge2Idx == SIZE_MAX)
443 continue;
444
445 Vec3D<NumericType> int1{}, int2{};
446 NumericType t1 = 0, t2 = 0;
447 bool cross1 = getSegmentIntersection(
448 pos0, pos1, mesh->nodes[mesh->lines[edge1Idx][0]],
449 mesh->nodes[mesh->lines[edge1Idx][1]], int1, t1);
450 bool cross2 = getSegmentIntersection(
451 pos0, pos1, mesh->nodes[mesh->lines[edge2Idx][0]],
452 mesh->nodes[mesh->lines[edge2Idx][1]], int2, t2);
453 if (!cross1 || !cross2)
454 continue;
455
456 if (t1 > t2) {
457 std::swap(t1, t2);
458 std::swap(int1, int2);
459 std::swap(edge1Idx, edge2Idx);
460 }
461
462 unsigned intNode1 = this->insertNode(int1);
463 unsigned intNode2 = this->insertNode(int2);
464
465 NumericType matM1 = currentMaterials[edge1Idx];
466 Vec3D<NumericType> normalM1 = currentNormals[edge1Idx];
467 NumericType matM2 = currentMaterials[edge2Idx];
468 Vec3D<NumericType> normalM2 = currentNormals[edge2Idx];
469 NumericType matL = currentMaterials[lineIdx];
470 Vec3D<NumericType> normalL = currentNormals[lineIdx];
471 unsigned e1n0 = mesh->lines[edge1Idx][0],
472 e1n1 = mesh->lines[edge1Idx][1];
473 unsigned e2n0 = mesh->lines[edge2Idx][0],
474 e2n1 = mesh->lines[edge2Idx][1];
475
476 uniqueElements.erase({(int)e1n0, (int)e1n1, 0});
477 mesh->lines[edge1Idx][1] = intNode1;
478 uniqueElements.insert({(int)e1n0, (int)intNode1, 0});
479 mesh->lines.push_back({intNode1, e1n1});
480 currentMaterials.push_back(matM1);
481 currentNormals.push_back(normalM1);
482 uniqueElements.insert({(int)intNode1, (int)e1n1, 0});
483
484 uniqueElements.erase({(int)e2n0, (int)e2n1, 0});
485 mesh->lines[edge2Idx][1] = intNode2;
486 uniqueElements.insert({(int)e2n0, (int)intNode2, 0});
487 mesh->lines.push_back({intNode2, e2n1});
488 currentMaterials.push_back(matM2);
489 currentNormals.push_back(normalM2);
490 uniqueElements.insert({(int)intNode2, (int)e2n1, 0});
491
492 uniqueElements.erase({(int)node0, (int)node1, 0});
493 mesh->lines[lineIdx][1] = intNode1;
494 uniqueElements.insert({(int)node0, (int)intNode1, 0});
495 mesh->lines.push_back({intNode2, node1});
496 currentMaterials.push_back(matL);
497 currentNormals.push_back(normalL);
498 uniqueElements.insert({(int)intNode2, (int)node1, 0});
499
500 return;
501 }
502 }
503 }
504
505 // Handle snapping/splitting of sharp corners between materials.
506 // activeGridIdx and lsLAtActive identify the active grid corner of the
507 // current cell (see findNearestLowerCorner) for the layer thickness check.
508 // cellIndices is the integer grid index of the current cell, used to restrict
509 // the double-crossing check to lower material corners in the same cell.
510 void
511 snapSharpCorners(unsigned l, size_t meshSizeBefore,
512 std::unordered_map<unsigned, std::vector<SharpCornerNode>>
513 &sharpCornerNodes,
514 const hrleIndex &activeGridIdx, NumericType lsLAtActive,
515 const hrleIndex &cellIndices) {
516 for (const auto &cornerPair : matSharpCornerNodes) {
517 unsigned cornerNodeId = cornerPair.first;
518 Vec3D<NumericType> cornerPos = cornerPair.second;
519
520 int snapToNodeId =
521 (l > 0) ? findNearestLowerCorner(l, cornerPos, sharpCornerNodes,
522 activeGridIdx, lsLAtActive)
523 : -1;
524
525 if (snapToNodeId >= 0) {
526 // Snap: replace cornerNodeId with snapToNodeId in recent elements
527 if (snapToNodeId != (int)cornerNodeId) {
528 if constexpr (D == 2) {
529 using I3 = typename ToSurfaceMesh<NumericType, D>::I3;
530 for (size_t i = mesh->lines.size(); i-- > meshSizeBefore;) {
531 for (int j = 0; j < 2; ++j) {
532 if (mesh->lines[i][j] == cornerNodeId)
533 mesh->lines[i][j] = snapToNodeId;
534 }
535
536 unsigned n0 = mesh->lines[i][0];
537 unsigned n1 = mesh->lines[i][1];
538 bool isDuplicate = uniqueElements.count({(int)n0, (int)n1, 0}) ||
539 uniqueElements.count({(int)n1, (int)n0, 0});
540
541 if (n0 == n1 || isDuplicate) {
542 mesh->lines.erase(mesh->lines.begin() + i);
543 currentMaterials.erase(currentMaterials.begin() + i);
544 currentNormals.erase(currentNormals.begin() + i);
545 } else {
546 uniqueElements.insert({(int)n0, (int)n1, 0});
547 }
548 }
549 }
550 }
551 sharpCornerNodes[l].push_back(
552 {static_cast<unsigned>(snapToNodeId), cornerPos});
553 } else {
554 // No snap - keep new corner and split any overlapping edge
555 sharpCornerNodes[l].push_back({cornerNodeId, cornerPos});
556 splitPreviousMaterialEdge(l, cornerNodeId, cornerPos, meshSizeBefore);
557 // For thin-layer case: the endpoints (nX, nY) of the polymer's sharp
558 // corner edges may lie on previous material surface edges.
559 // Collect all endpoints first, then apply all splits before crossing
560 // checks. Interleaving the two would risk index-shifting: erases inside
561 // splitPreviousMaterialEdge's duplicate-removal loop can invalidate
562 // pIdx before handleCellCornerCrossings is called.
563 size_t sizeNow = mesh->lines.size();
564 std::vector<std::pair<unsigned, Vec3D<NumericType>>> endpointsToSplit;
565 for (size_t pIdx = meshSizeBefore; pIdx < sizeNow; ++pIdx) {
566 if (static_cast<unsigned>(currentMaterials[pIdx]) != l)
567 continue;
568 unsigned n0 = mesh->lines[pIdx][0];
569 unsigned n1 = mesh->lines[pIdx][1];
570 if (n0 != cornerNodeId)
571 endpointsToSplit.push_back({n0, mesh->nodes[n0]});
572 if (n1 != cornerNodeId)
573 endpointsToSplit.push_back({n1, mesh->nodes[n1]});
574 }
575 for (auto &[nodeId, pos] : endpointsToSplit)
576 splitPreviousMaterialEdge(l, nodeId, pos, meshSizeBefore);
577 // Run crossing checks after all lower-material splits are complete.
578 sizeNow = mesh->lines.size();
579 for (size_t pIdx = meshSizeBefore; pIdx < sizeNow; ++pIdx) {
580 if (static_cast<unsigned>(currentMaterials[pIdx]) != l)
581 continue;
582 handleCellCornerCrossings(pIdx, sharpCornerNodes, l, cellIndices);
583 }
584 }
585 }
586 }
587
588public:
589 ToMultiSurfaceMesh(double minNodeDistFactor = 0.01, double eps = 1e-12)
590 : ToSurfaceMesh<NumericType, D>(minNodeDistFactor, eps) {}
591
592 ToMultiSurfaceMesh(SmartPointer<lsDomainType> passedLevelSet,
593 SmartPointer<viennals::Mesh<NumericType>> passedMesh,
594 double minNodeDistFactor = 0.01, double eps = 1e-12)
595 : ToSurfaceMesh<NumericType, D>(passedLevelSet, passedMesh,
596 minNodeDistFactor, eps) {}
597
599 std::vector<SmartPointer<lsDomainType>> const &passedLevelSets,
600 SmartPointer<viennals::Mesh<NumericType>> passedMesh,
601 double minNodeDistFactor = 0.01, double eps = 1e-12)
602 : ToSurfaceMesh<NumericType, D>(minNodeDistFactor, eps) {
603 levelSets = passedLevelSets;
604 mesh = passedMesh;
605 }
606
608 double minNodeDistFactor = 0.01, double eps = 1e-12)
609 : ToSurfaceMesh<NumericType, D>(minNodeDistFactor, eps) {
610 mesh = passedMesh;
611 }
612
613 void insertNextLevelSet(SmartPointer<lsDomainType> passedLevelSet) {
614 levelSets.push_back(passedLevelSet);
615 }
616
617 void clearLevelSets() { levelSets.clear(); }
618
619 void setMaterialMap(SmartPointer<MaterialMap> passedMaterialMap) {
620 materialMap = passedMaterialMap;
621 }
622
623 void apply() override {
624 if (levelSets.empty()) {
625 Logger::getInstance()
626 .addError("No level set was passed to ToMultiSurfaceMesh.")
627 .print();
628 return;
629 }
630 if (mesh == nullptr) {
631 Logger::getInstance()
632 .addError("No mesh was passed to ToMultiSurfaceMesh.")
633 .print();
634 return;
635 }
636
637 mesh->clear();
638 currentGridDelta = levelSets.front()->getGrid().getGridDelta();
639 for (unsigned i = 0; i < D; ++i) {
640 mesh->minimumExtent[i] = std::numeric_limits<NumericType>::max();
641 mesh->maximumExtent[i] = std::numeric_limits<NumericType>::lowest();
642 }
643
644 using nodeContainerType = std::map<viennahrle::Index<D>, unsigned>;
645
646 nodeContainerType nodes[D];
647 nodeContainerType faceNodes[D];
648 nodeContainerType cornerNodes;
649
650 nodeIdByBin.clear();
651 uniqueElements.clear();
652 currentNormals.clear();
653 currentMaterials.clear();
654
655 const bool useMaterialMap = materialMap != nullptr;
656 const bool sharpCorners = generateSharpCorners;
657
658 if (sharpCorners && D == 3) {
659 Logger::getInstance()
660 .addWarning("Sharp corner generation in 3D is experimental and may "
661 "produce suboptimal meshes. Use with caution.")
662 .print();
663 }
664
665 // an iterator for each level set
666 std::vector<viennahrle::ConstSparseCellIterator<hrleDomainType>> cellIts;
667 for (const auto &ls : levelSets)
668 cellIts.emplace_back(ls->getDomain());
669
670 // Explicit storage for sharp corner nodes, keyed by material index
671 std::unordered_map<unsigned, std::vector<SharpCornerNode>> sharpCornerNodes;
672
673 for (unsigned l = 0; l < levelSets.size(); l++) {
675 if (useMaterialMap) {
676 currentMaterialId = materialMap->getMaterialId(l);
677 } else {
678 currentMaterialId = static_cast<unsigned>(l);
679 }
680
681 normalVectorData = nullptr;
682 if (sharpCorners) {
685 normalCalculator.setMethod(
687 normalCalculator.setMaxValue(std::numeric_limits<NumericType>::max());
688 normalCalculator.apply();
689 normalVectorData = currentLevelSet->getPointData().getVectorData(
691 }
692
693 viennahrle::ConstSparseIterator<hrleDomainType> valueIt(
694 currentLevelSet->getDomain());
695
696 // iterate over all active surface points
697 for (auto cellIt = cellIts[l]; !cellIt.isFinished(); cellIt.next()) {
698 for (int u = 0; u < D; u++) {
699 while (!nodes[u].empty() &&
700 nodes[u].begin()->first <
701 viennahrle::Index<D>(cellIt.getIndices()))
702 nodes[u].erase(nodes[u].begin());
703
704 if constexpr (D == 3) {
705 while (!faceNodes[u].empty() &&
706 faceNodes[u].begin()->first <
707 viennahrle::Index<D>(cellIt.getIndices()))
708 faceNodes[u].erase(faceNodes[u].begin());
709
710 if (u == 0) {
711 while (!cornerNodes.empty() &&
712 cornerNodes.begin()->first <
713 viennahrle::Index<D>(cellIt.getIndices()))
714 cornerNodes.erase(cornerNodes.begin());
715 }
716 }
717 }
718
719 unsigned signs = 0;
720 for (int i = 0; i < (1 << D); i++) {
721 if (cellIt.getCorner(i).getValue() >= NumericType(0))
722 signs |= (1 << i);
723 }
724
725 // all corners have the same sign, so no surface here
726 if (signs == 0)
727 continue;
728 if (signs == (1 << (1 << D)) - 1)
729 continue;
730
731 // Check if this cell is at a material boundary with the immediate
732 // material below
733 bool atMaterialBoundary = false;
734 unsigned touchingMaterial = 0; // Which material we're touching
735 if (sharpCorners && l > 0) {
736 touchingMaterial = l - 1;
737 viennahrle::ConstSparseIterator<hrleDomainType> prevIt(
738 levelSets[touchingMaterial]->getDomain());
739 for (int i = 0; i < (1 << D); ++i) {
740 hrleIndex cornerIdx =
741 cellIt.getIndices() + viennahrle::BitMaskToIndex<D>(i);
742 prevIt.goToIndices(cornerIdx);
743 if (prevIt.getValue() <=
744 NumericType(0)) { // Inside previous material
745 atMaterialBoundary = true;
746 break;
747 }
748 }
749 }
750
751 // Attempt to generate sharp features if enabled
752 bool perfectCornerFound = false;
753 if (sharpCorners) {
754 int countNeg = 0;
755 int countPos = 0;
756 int negMask = 0;
757 int posMask = 0;
758 for (int i = 0; i < (1 << D); ++i) {
759 NumericType val = cellIt.getCorner(i).getValue();
760 if (val < -epsilon) {
761 countNeg++;
762 negMask |= (1 << i);
763 } else if (val > epsilon) {
764 countPos++;
765 posMask |= (1 << i);
766 } else {
767 if (val >= NumericType(0)) {
768 countPos++;
769 posMask |= (1 << i);
770 } else {
771 countNeg++;
772 negMask |= (1 << i);
773 }
774 }
775 }
776
777 // Clear the base class's recent corner nodes before generation
778 matSharpCornerNodes.clear();
779
780 // Track mesh size before generating sharp corners
781 size_t meshSizeBefore = 0;
782 if constexpr (D == 2) {
783 meshSizeBefore = mesh->lines.size();
784 perfectCornerFound =
785 this->generateSharpCorner2D(cellIt, countNeg, countPos, negMask,
786 posMask, nodes, nullptr, valueIt);
787 } else if constexpr (D == 3) {
788 if (countNeg == 2 || countPos == 2) {
789 perfectCornerFound = this->generateSharpEdge3D(
790 cellIt, countNeg, countPos, negMask, posMask, nodes,
791 faceNodes, nullptr, valueIt);
792 } else if (countNeg == 1 || countPos == 1) {
793 perfectCornerFound = this->generateSharpCorner3D(
794 cellIt, countNeg, countPos, negMask, posMask, nodes,
795 cornerNodes, faceNodes, nullptr, valueIt);
796 } else if (countNeg == 3 || countPos == 3) {
797 perfectCornerFound = this->generateSharpL3D(
798 cellIt, countNeg, countPos, negMask, posMask, nodes,
799 faceNodes, nullptr, valueIt);
800 }
801 }
802
803 // If sharp corners were generated, check if they should snap to
804 // existing corners or split edges of previous materials.
805 // Identify the "active" (single-minority sign) grid corner so the
806 // layer thickness check inside snapSharpCorners can use grid-level LS
807 // values instead of unreliable sub-grid interpolation.
808 if (D == 2 && perfectCornerFound && !matSharpCornerNodes.empty()) {
809 int activeBit = -1;
810 if (countNeg == 1) {
811 for (int i = 0; i < (1 << D); ++i)
812 if ((negMask >> i) & 1) {
813 activeBit = i;
814 break;
815 }
816 } else if (countPos == 1) {
817 for (int i = 0; i < (1 << D); ++i)
818 if ((posMask >> i) & 1) {
819 activeBit = i;
820 break;
821 }
822 }
823 hrleIndex activeGridIdx = cellIt.getIndices();
824 NumericType lsLAtActive = NumericType(0);
825 if (activeBit >= 0) {
826 activeGridIdx += viennahrle::BitMaskToIndex<D>(activeBit);
827 lsLAtActive = cellIt.getCorner(activeBit).getValue();
828 }
829 snapSharpCorners(l, meshSizeBefore, sharpCornerNodes, activeGridIdx,
830 lsLAtActive, cellIt.getIndices());
831 }
832 }
833
834 if (perfectCornerFound)
835 continue;
836
837 if constexpr (D == 3) {
838 // Stitch to perfect corners/edges
839 for (int axis = 0; axis < 3; ++axis) {
840 for (int d = 0; d < 2; ++d) {
841 viennahrle::Index<D> faceKey(cellIt.getIndices());
842 if (d == 1)
843 faceKey[axis]++;
844
845 auto it = faceNodes[axis].find(faceKey);
846 if (it != faceNodes[axis].end()) {
847 const unsigned faceNodeId = it->second;
848 const int *Triangles =
850
851 for (; Triangles[0] != -1; Triangles += 3) {
852 std::vector<unsigned> face_edge_nodes;
853 for (int i = 0; i < 3; ++i) {
854 int edge = Triangles[i];
855 int c0 = this->corner0[edge];
856 int c1 = this->corner1[edge];
857 bool onFace =
858 (((c0 >> axis) & 1) == d) && (((c1 >> axis) & 1) == d);
859 if (onFace) {
860 face_edge_nodes.push_back(
861 this->getNode(cellIt, edge, nodes, nullptr));
862 }
863 }
864 if (face_edge_nodes.size() == 2) {
865 this->insertElement(std::array<unsigned, 3>{
866 face_edge_nodes[0], face_edge_nodes[1], faceNodeId});
867 }
868 }
869 }
870 }
871 }
872 }
873
874 // for each element
875 const int *Triangles;
876 if constexpr (D == 2) {
878 } else {
880 }
881
882 for (; Triangles[0] != -1; Triangles += D) {
883 std::array<unsigned, D> nodeNumbers;
884
885 // for each node
886 for (int n = 0; n < D; n++) {
887 const int edge = Triangles[n];
888
889 unsigned p0 = this->corner0[edge];
890
891 // determine direction of edge
892 unsigned dir = this->direction[edge];
893
894 // look for existing surface node
895 viennahrle::Index<D> d(cellIt.getIndices());
896 auto p0B = viennahrle::BitMaskToIndex<D>(p0);
897 d += p0B;
898
899 auto nodeIt = nodes[dir].find(d);
900 if (nodeIt != nodes[dir].end()) {
901 nodeNumbers[n] = nodeIt->second;
902
903 // Check if this cached node is a sharp corner that should be
904 // inherited
905 if (sharpCorners && l > 0 && atMaterialBoundary) {
906 unsigned cachedNodeId = nodeIt->second;
907 // Check if this is a corner from the material below
908 auto it = sharpCornerNodes.find(touchingMaterial);
909 if (it != sharpCornerNodes.end()) {
910 for (const auto &scn : it->second) {
911 if (scn.nodeId == cachedNodeId) {
912 // Only inherit if not already inherited
913 auto &lCorners = sharpCornerNodes[l];
914 bool alreadyInherited = false;
915 for (const auto &existing : lCorners) {
916 if (existing.nodeId == cachedNodeId) {
917 alreadyInherited = true;
918 break;
919 }
920 }
921 if (!alreadyInherited)
922 lCorners.push_back({cachedNodeId, scn.position});
923 break;
924 }
925 }
926 }
927 }
928 } else {
929 // Node does not exist yet - try snapping to a sharp corner
930 // from the touching material at material boundaries
931 int snapNodeId = -1;
932 Vec3D<NumericType> snapPos{};
933 if (sharpCorners && l > 0 && atMaterialBoundary) {
934 auto tmIt = sharpCornerNodes.find(touchingMaterial);
935 if (tmIt != sharpCornerNodes.end() && !tmIt->second.empty()) {
936 Vec3D<NumericType> nodePos =
937 this->computeNodePosition(cellIt, edge);
938
939 NumericType minDist2 =
941 for (const auto &scn : tmIt->second) {
942 NumericType dist2 = NumericType(0);
943 for (int i = 0; i < D; ++i) {
944 NumericType dd = scn.position[i] - nodePos[i];
945 dist2 += dd * dd;
946 }
947 if (dist2 < minDist2) {
948 minDist2 = dist2;
949 snapNodeId = scn.nodeId;
950 snapPos = scn.position;
951 }
952 }
953 }
954 }
955
956 if (snapNodeId >= 0) {
957 nodeNumbers[n] = snapNodeId;
958 nodes[dir][d] = snapNodeId; // Cache it
959 // Inherit this corner so upper materials can find it
960 sharpCornerNodes[l].push_back(
961 {static_cast<unsigned>(snapNodeId), snapPos});
962 } else {
963 nodeNumbers[n] = this->getNode(cellIt, edge, nodes, nullptr);
964 }
965 }
966 }
967
968 // Check for duplicates before inserting (standard Marching Cubes)
969 bool isDuplicate = false;
970 if constexpr (D == 2) {
971 using I3 = typename ToSurfaceMesh<NumericType, D>::I3;
972 if (uniqueElements.find({(int)nodeNumbers[1], (int)nodeNumbers[0],
973 0}) != uniqueElements.end())
974 isDuplicate = true;
975 }
976
977 if (!isDuplicate) {
978 this->insertElement(nodeNumbers);
979 if (sharpCorners && l > 0)
980 handleEdgeCornerInteractions(mesh->lines.size() - 1,
981 sharpCornerNodes, l,
982 cellIt.getIndices());
983 }
984 }
985 }
986 }
987
988 this->scaleMesh();
989
990 mesh->cellData.insertNextScalarData(currentMaterials, "MaterialIds");
991 mesh->cellData.insertNextVectorData(currentNormals, "Normals");
992 mesh->triangles.shrink_to_fit();
993 mesh->nodes.shrink_to_fit();
994 }
995};
996
997PRECOMPILE_PRECISION_DIMENSION(ToMultiSurfaceMesh)
998
999} // namespace viennals
float NumericType
Definition AirGapDeposition.cpp:24
constexpr int D
Definition Epitaxy.cpp:11
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< NumericType, D > DomainType
Definition lsDomain.hpp:33
This class holds an explicit mesh, which is always given in 3 dimensions. If it describes a 2D mesh,...
Definition lsMesh.hpp:22
ToMultiSurfaceMesh(SmartPointer< lsDomainType > passedLevelSet, SmartPointer< viennals::Mesh< NumericType > > passedMesh, double minNodeDistFactor=0.01, double eps=1e-12)
Definition lsToMultiSurfaceMesh.hpp:592
ToMultiSurfaceMesh(std::vector< SmartPointer< lsDomainType > > const &passedLevelSets, SmartPointer< viennals::Mesh< NumericType > > passedMesh, double minNodeDistFactor=0.01, double eps=1e-12)
Definition lsToMultiSurfaceMesh.hpp:598
void apply() override
Definition lsToMultiSurfaceMesh.hpp:623
void setMaterialMap(SmartPointer< MaterialMap > passedMaterialMap)
Definition lsToMultiSurfaceMesh.hpp:619
ToMultiSurfaceMesh(double minNodeDistFactor=0.01, double eps=1e-12)
Definition lsToMultiSurfaceMesh.hpp:589
void clearLevelSets()
Definition lsToMultiSurfaceMesh.hpp:617
void insertNextLevelSet(SmartPointer< lsDomainType > passedLevelSet)
Definition lsToMultiSurfaceMesh.hpp:613
ToMultiSurfaceMesh(SmartPointer< viennals::Mesh< NumericType > > passedMesh, double minNodeDistFactor=0.01, double eps=1e-12)
Definition lsToMultiSurfaceMesh.hpp:607
SmartPointer< Mesh< NumericType > > mesh
Definition lsToSurfaceMesh.hpp:33
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
static constexpr unsigned int direction[12]
Definition lsToSurfaceMesh.hpp:77
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::unordered_map< I3, unsigned, I3Hash > nodeIdByBin
Definition lsToSurfaceMesh.hpp:61
const NumericType epsilon
Definition lsToSurfaceMesh.hpp:35
Vec3D< NumericType > 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::vector< std::pair< unsigned, Vec3D< NumericType > > > matSharpCornerNodes
Definition lsToSurfaceMesh.hpp:71
static constexpr unsigned int corner1[12]
Definition lsToSurfaceMesh.hpp:75
const double minNodeDistanceFactor
Definition lsToSurfaceMesh.hpp:36
void insertElement(const std::array< unsigned, D > &nodeNumbers)
Definition lsToSurfaceMesh.hpp:1654
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< NumericType > const &pos)
Definition lsToSurfaceMesh.hpp:1692
PointData< NumericType >::VectorDataType * normalVectorData
Definition lsToSurfaceMesh.hpp:68
std::vector< Vec3D< NumericType > > currentNormals
Definition lsToSurfaceMesh.hpp:63
void scaleMesh()
Definition lsToSurfaceMesh.hpp:342
std::vector< NumericType > currentMaterials
Definition lsToSurfaceMesh.hpp:64
NumericType currentMaterialId
Definition lsToSurfaceMesh.hpp:66
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:40