ViennaLS
Loading...
Searching...
No Matches
lsToMultiSurfaceMesh.hpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <iostream>
6#include <map>
7#include <unordered_map>
8#include <unordered_set>
9
10#include <hrleSparseCellIterator.hpp>
11#include <lsDomain.hpp>
12#include <lsExpand.hpp>
13#include <lsMarchingCubes.hpp>
14#include <lsMaterialMap.hpp>
15#include <lsMesh.hpp>
16
17namespace viennals {
18
19using namespace viennacore;
20
21template <class NumericType, int D, class T = NumericType>
23
24 typedef viennals::Domain<NumericType, D> lsDomainType;
25 typedef typename viennals::Domain<NumericType, D>::DomainType hrleDomainType;
26
27 std::vector<SmartPointer<lsDomainType>> levelSets;
28 SmartPointer<viennals::Mesh<T>> mesh = nullptr;
29 SmartPointer<MaterialMap> materialMap = nullptr;
30
31 const double epsilon;
32 const double minNodeDistanceFactor;
33
34 struct I3 {
35 int x, y, z;
36 bool operator==(const I3 &o) const {
37 return x == o.x && y == o.y && z == o.z;
38 }
39 };
40
41 struct I3Hash {
42 size_t operator()(const I3 &k) const {
43 // 64-bit mix
44 uint64_t a = (uint64_t)(uint32_t)k.x;
45 uint64_t b = (uint64_t)(uint32_t)k.y;
46 uint64_t c = (uint64_t)(uint32_t)k.z;
47 uint64_t h = a * 0x9E3779B185EBCA87ULL;
48 h ^= b + 0xC2B2AE3D27D4EB4FULL + (h << 6) + (h >> 2);
49 h ^= c + 0x165667B19E3779F9ULL + (h << 6) + (h >> 2);
50 return (size_t)h;
51 }
52 };
53
54public:
55 ToMultiSurfaceMesh(double eps = 1e-12, double minNodeDistFactor = 0.05)
56 : epsilon(eps), minNodeDistanceFactor(minNodeDistFactor) {}
57
58 ToMultiSurfaceMesh(SmartPointer<lsDomainType> passedLevelSet,
59 SmartPointer<viennals::Mesh<T>> passedMesh,
60 double eps = 1e-12, double minNodeDistFactor = 0.05)
61 : mesh(passedMesh), epsilon(eps),
62 minNodeDistanceFactor(minNodeDistFactor) {
63 levelSets.push_back(passedLevelSet);
64 }
65
66 ToMultiSurfaceMesh(SmartPointer<viennals::Mesh<T>> passedMesh,
67 double eps = 1e-12, double minNodeDistFactor = 0.05)
68 : mesh(passedMesh), epsilon(eps),
69 minNodeDistanceFactor(minNodeDistFactor) {}
70
71 void setMesh(SmartPointer<viennals::Mesh<T>> passedMesh) {
72 mesh = passedMesh;
73 }
74
75 void insertNextLevelSet(SmartPointer<lsDomainType> passedLevelSet) {
76 levelSets.push_back(passedLevelSet);
77 }
78
79 void clearLevelSets() { levelSets.clear(); }
80
81 void setMaterialMap(SmartPointer<MaterialMap> passedMaterialMap) {
82 materialMap = passedMaterialMap;
83 }
84
85 void apply() {
86 if (levelSets.empty()) {
87 Logger::getInstance()
88 .addError("No level set was passed to ToMultiSurfaceMesh.")
89 .print();
90 return;
91 }
92 if (mesh == nullptr) {
93 Logger::getInstance()
94 .addError("No mesh was passed to ToMultiSurfaceMesh.")
95 .print();
96 return;
97 }
98
99 mesh->clear();
100 const auto gridDelta = levelSets.front()->getGrid().getGridDelta();
101 for (unsigned i = 0; i < D; ++i) {
102 mesh->minimumExtent[i] = std::numeric_limits<NumericType>::max();
103 mesh->maximumExtent[i] = std::numeric_limits<NumericType>::lowest();
104 }
105
106 constexpr unsigned int corner0[12] = {0, 1, 2, 0, 4, 5, 6, 4, 0, 1, 3, 2};
107 constexpr unsigned int corner1[12] = {1, 3, 3, 2, 5, 7, 7, 6, 4, 5, 7, 6};
108 constexpr unsigned int direction[12] = {0, 1, 0, 1, 0, 1, 0, 1, 2, 2, 2, 2};
109
110 typedef std::map<viennahrle::Index<D>, unsigned> nodeContainerType;
111
112 nodeContainerType nodes[D];
113 const double minNodeDistance = gridDelta * minNodeDistanceFactor;
114 std::unordered_map<I3, unsigned, I3Hash> nodeIdByBin;
115 std::unordered_set<I3, I3Hash> uniqueElements;
116
117 typename nodeContainerType::iterator nodeIt;
118
119 std::vector<Vec3D<T>> normals;
120 std::vector<T> materials;
121
122 const bool checkNodeFlag = minNodeDistanceFactor > 0;
123 const bool useMaterialMap = materialMap != nullptr;
124
125 auto quantize = [&](const Vec3D<T> &p) -> I3 {
126 const T inv = T(1) / minNodeDistance;
127 return {(int)std::llround(p[0] * inv), (int)std::llround(p[1] * inv),
128 (int)std::llround(p[2] * inv)};
129 };
130
131 auto elementI3 = [&](const std::array<unsigned, D> &element) -> I3 {
132 if constexpr (D == 2)
133 return {(int)element[0], (int)element[1], 0};
134 else
135 return {(int)element[0], (int)element[1], (int)element[2]};
136 };
137
138 // an iterator for each level set
139 std::vector<viennahrle::ConstSparseCellIterator<hrleDomainType>> cellIts;
140 for (const auto &ls : levelSets)
141 cellIts.emplace_back(ls->getDomain());
142
143 for (unsigned l = 0; l < levelSets.size(); l++) {
144 // iterate over all active surface points
145 for (auto cellIt = cellIts[l]; !cellIt.isFinished(); cellIt.next()) {
146 for (int u = 0; u < D; u++) {
147 while (!nodes[u].empty() &&
148 nodes[u].begin()->first <
149 viennahrle::Index<D>(cellIt.getIndices()))
150 nodes[u].erase(nodes[u].begin());
151 }
152
153 unsigned signs = 0;
154 for (int i = 0; i < (1 << D); i++) {
155 if (cellIt.getCorner(i).getValue() >= NumericType(0))
156 signs |= (1 << i);
157 }
158
159 // all corners have the same sign, so no surface here
160 if (signs == 0)
161 continue;
162 if (signs == (1 << (1 << D)) - 1)
163 continue;
164
165 // for each element
166 const int *Triangles;
167 if constexpr (D == 2) {
169 } else {
171 }
172
173 for (; Triangles[0] != -1; Triangles += D) {
174 std::array<unsigned, D> nod_numbers;
175
176 // for each node
177 for (int n = 0; n < D; n++) {
178 const int edge = Triangles[n];
179
180 unsigned p0 = corner0[edge];
181 unsigned p1 = corner1[edge];
182
183 // determine direction of edge
184 unsigned dir = direction[edge];
185
186 // look for existing surface node
187 viennahrle::Index<D> d(cellIt.getIndices());
188 auto p0B = viennahrle::BitMaskToIndex<D>(p0);
189 d += p0B;
190
191 nodeIt = nodes[dir].find(d);
192 if (nodeIt != nodes[dir].end()) {
193 nod_numbers[n] = nodeIt->second;
194 } else {
195 // if node does not exist yet
196 // calculate coordinate of new node
197 Vec3D<T> cc{}; // initialise with zeros
198 for (int z = 0; z < D; z++) {
199 if (z != dir) {
200 // TODO might not need BitMaskToVector here, just check if z
201 // bit is set
202 cc[z] = static_cast<T>(cellIt.getIndices(z) + p0B[z]);
203 } else {
204 auto d0 = static_cast<T>(cellIt.getCorner(p0).getValue());
205 auto d1 = static_cast<T>(cellIt.getCorner(p1).getValue());
206
207 // calculate the surface-grid intersection point
208 if (d0 == -d1) { // includes case where d0=d1=0
209 cc[z] = static_cast<T>(cellIt.getIndices(z)) + 0.5;
210 } else {
211 if (std::abs(d0) <= std::abs(d1)) {
212 cc[z] = static_cast<T>(cellIt.getIndices(z)) +
213 (d0 / (d0 - d1));
214 } else {
215 cc[z] = static_cast<T>(cellIt.getIndices(z) + 1) -
216 (d1 / (d1 - d0));
217 }
218 }
219 cc[z] = std::max(
220 cc[z], static_cast<T>(cellIt.getIndices(z) + epsilon));
221 cc[z] = std::min(
222 cc[z],
223 static_cast<T>((cellIt.getIndices(z) + 1) - epsilon));
224 }
225 cc[z] *= gridDelta;
226 }
227
228 int nodeIdx = -1;
229 if (checkNodeFlag) {
230 auto q = quantize(cc);
231 auto it = nodeIdByBin.find(q);
232 if (it != nodeIdByBin.end())
233 nodeIdx = it->second;
234 }
235 if (nodeIdx >= 0) {
236 nod_numbers[n] = nodeIdx;
237 } else {
238 // insert new node
239 nod_numbers[n] =
240 mesh->insertNextNode(cc); // insert new surface node
241 nodes[dir][d] = nod_numbers[n];
242 if (checkNodeFlag)
243 nodeIdByBin.emplace(quantize(cc), nod_numbers[n]);
244
245 for (int a = 0; a < D; a++) {
246 if (cc[a] < mesh->minimumExtent[a])
247 mesh->minimumExtent[a] = cc[a];
248 if (cc[a] > mesh->maximumExtent[a])
249 mesh->maximumExtent[a] = cc[a];
250 }
251 }
252 }
253 }
254
255 // check for degenerate triangles
256 if (!triangleMisformed(nod_numbers)) {
257
258 // only insert if not already present
259 if (uniqueElements.insert(elementI3(nod_numbers)).second) {
260 Vec3D<T> normal;
261 if constexpr (D == 2) {
262 normal = Vec3D<T>{-(mesh->nodes[nod_numbers[1]][1] -
263 mesh->nodes[nod_numbers[0]][1]),
264 mesh->nodes[nod_numbers[1]][0] -
265 mesh->nodes[nod_numbers[0]][0],
266 T(0)};
267 } else {
268 normal = calculateNormal(mesh->nodes[nod_numbers[0]],
269 mesh->nodes[nod_numbers[1]],
270 mesh->nodes[nod_numbers[2]]);
271 }
272
273 double n2 = normal[0] * normal[0] + normal[1] * normal[1] +
274 normal[2] * normal[2];
275 if (n2 > epsilon) {
276 // insert new element
277 mesh->insertNextElement(nod_numbers);
278
279 // insert normal
280 T invn = static_cast<T>(1.) / std::sqrt(static_cast<T>(n2));
281 for (int d = 0; d < D; d++) {
282 normal[d] *= invn;
283 }
284 normals.push_back(normal);
285
286 // insert material
287 if (useMaterialMap) {
288 materials.push_back(materialMap->getMaterialId(l));
289 } else {
290 materials.push_back(static_cast<T>(l));
291 }
292 }
293 }
294 }
295 }
296 }
297 }
298
299 mesh->cellData.insertNextScalarData(materials, "MaterialIds");
300 mesh->cellData.insertNextVectorData(normals, "Normals");
301 mesh->triangles.shrink_to_fit();
302 mesh->nodes.shrink_to_fit();
303 }
304
305private:
306 static inline bool
307 triangleMisformed(const std::array<unsigned, D> &nod_numbers) noexcept {
308 if constexpr (D == 3) {
309 return nod_numbers[0] == nod_numbers[1] ||
310 nod_numbers[0] == nod_numbers[2] ||
311 nod_numbers[1] == nod_numbers[2];
312 } else {
313 return nod_numbers[0] == nod_numbers[1];
314 }
315 }
316
317 static inline Vec3D<NumericType>
318 calculateNormal(const Vec3D<NumericType> &nodeA,
319 const Vec3D<NumericType> &nodeB,
320 const Vec3D<NumericType> &nodeC) noexcept {
321 return CrossProduct(nodeB - nodeA, nodeC - nodeA);
322 }
323};
324
326
327} // namespace viennals
float NumericType
Definition AirGapDeposition.cpp:24
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
Class containing all information about the level set, including the dimensions of the domain,...
Definition lsDomain.hpp:27
viennahrle::Domain< T, D > DomainType
Definition lsDomain.hpp:32
This class holds an explicit mesh, which is always given in 3 dimensions. If it describes a 2D mesh,...
Definition lsMesh.hpp:22
Definition lsToMultiSurfaceMesh.hpp:22
void setMaterialMap(SmartPointer< MaterialMap > passedMaterialMap)
Definition lsToMultiSurfaceMesh.hpp:81
void setMesh(SmartPointer< viennals::Mesh< T > > passedMesh)
Definition lsToMultiSurfaceMesh.hpp:71
void apply()
Definition lsToMultiSurfaceMesh.hpp:85
void clearLevelSets()
Definition lsToMultiSurfaceMesh.hpp:79
ToMultiSurfaceMesh(double eps=1e-12, double minNodeDistFactor=0.05)
Definition lsToMultiSurfaceMesh.hpp:55
ToMultiSurfaceMesh(SmartPointer< lsDomainType > passedLevelSet, SmartPointer< viennals::Mesh< T > > passedMesh, double eps=1e-12, double minNodeDistFactor=0.05)
Definition lsToMultiSurfaceMesh.hpp:58
void insertNextLevelSet(SmartPointer< lsDomainType > passedLevelSet)
Definition lsToMultiSurfaceMesh.hpp:75
ToMultiSurfaceMesh(SmartPointer< viennals::Mesh< T > > passedMesh, double eps=1e-12, double minNodeDistFactor=0.05)
Definition lsToMultiSurfaceMesh.hpp:66
#define PRECOMPILE_PRECISION_DIMENSION(className)
Definition lsPreCompileMacros.hpp:24
Definition lsAdvect.hpp:36