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 setMaterialMap(SmartPointer<MaterialMap> passedMaterialMap) {
80 materialMap = passedMaterialMap;
81 }
82
83 void apply() {
84 if (levelSets.empty()) {
85 Logger::getInstance()
86 .addError("No level set was passed to ToMultiSurfaceMesh.")
87 .print();
88 return;
89 }
90 if (mesh == nullptr) {
91 Logger::getInstance()
92 .addError("No mesh was passed to ToMultiSurfaceMesh.")
93 .print();
94 return;
95 }
96
97 mesh->clear();
98 const auto gridDelta = levelSets.front()->getGrid().getGridDelta();
99 for (unsigned i = 0; i < D; ++i) {
100 mesh->minimumExtent[i] = std::numeric_limits<NumericType>::max();
101 mesh->maximumExtent[i] = std::numeric_limits<NumericType>::lowest();
102 }
103
104 constexpr unsigned int corner0[12] = {0, 1, 2, 0, 4, 5, 6, 4, 0, 1, 3, 2};
105 constexpr unsigned int corner1[12] = {1, 3, 3, 2, 5, 7, 7, 6, 4, 5, 7, 6};
106 constexpr unsigned int direction[12] = {0, 1, 0, 1, 0, 1, 0, 1, 2, 2, 2, 2};
107
108 typedef std::map<viennahrle::Index<D>, unsigned> nodeContainerType;
109
110 nodeContainerType nodes[D];
111 const double minNodeDistance = gridDelta * minNodeDistanceFactor;
112 std::unordered_map<I3, unsigned, I3Hash> nodeIdByBin;
113 std::unordered_set<I3, I3Hash> uniqueElements;
114
115 typename nodeContainerType::iterator nodeIt;
116
117 std::vector<Vec3D<T>> normals;
118 std::vector<T> materials;
119
120 const bool checkNodeFlag = minNodeDistanceFactor > 0;
121 const bool useMaterialMap = materialMap != nullptr;
122
123 auto quantize = [&](const Vec3D<T> &p) -> I3 {
124 const T inv = T(1) / minNodeDistance;
125 return {(int)std::llround(p[0] * inv), (int)std::llround(p[1] * inv),
126 (int)std::llround(p[2] * inv)};
127 };
128
129 auto elementI3 = [&](const std::array<unsigned, D> &element) -> I3 {
130 if constexpr (D == 2)
131 return {(int)element[0], (int)element[1], 0};
132 else
133 return {(int)element[0], (int)element[1], (int)element[2]};
134 };
135
136 // an iterator for each level set
137 std::vector<viennahrle::ConstSparseCellIterator<hrleDomainType>> cellIts;
138 for (const auto &ls : levelSets)
139 cellIts.emplace_back(ls->getDomain());
140
141 for (unsigned l = 0; l < levelSets.size(); l++) {
142 // iterate over all active surface points
143 for (auto cellIt = cellIts[l]; !cellIt.isFinished(); cellIt.next()) {
144 for (int u = 0; u < D; u++) {
145 while (!nodes[u].empty() &&
146 nodes[u].begin()->first <
147 viennahrle::Index<D>(cellIt.getIndices()))
148 nodes[u].erase(nodes[u].begin());
149 }
150
151 unsigned signs = 0;
152 for (int i = 0; i < (1 << D); i++) {
153 if (cellIt.getCorner(i).getValue() >= NumericType(0))
154 signs |= (1 << i);
155 }
156
157 // all corners have the same sign, so no surface here
158 if (signs == 0)
159 continue;
160 if (signs == (1 << (1 << D)) - 1)
161 continue;
162
163 // for each element
164 const int *Triangles;
165 if constexpr (D == 2) {
167 } else {
169 }
170
171 for (; Triangles[0] != -1; Triangles += D) {
172 std::array<unsigned, D> nod_numbers;
173
174 // for each node
175 for (int n = 0; n < D; n++) {
176 const int edge = Triangles[n];
177
178 unsigned p0 = corner0[edge];
179 unsigned p1 = corner1[edge];
180
181 // determine direction of edge
182 unsigned dir = direction[edge];
183
184 // look for existing surface node
185 viennahrle::Index<D> d(cellIt.getIndices());
186 auto p0B = viennahrle::BitMaskToIndex<D>(p0);
187 d += p0B;
188
189 nodeIt = nodes[dir].find(d);
190 if (nodeIt != nodes[dir].end()) {
191 nod_numbers[n] = nodeIt->second;
192 } else {
193 // if node does not exist yet
194 // calculate coordinate of new node
195 Vec3D<T> cc{}; // initialise with zeros
196 for (int z = 0; z < D; z++) {
197 if (z != dir) {
198 // TODO might not need BitMaskToVector here, just check if z
199 // bit is set
200 cc[z] = static_cast<T>(cellIt.getIndices(z) + p0B[z]);
201 } else {
202 auto d0 = static_cast<T>(cellIt.getCorner(p0).getValue());
203 auto d1 = static_cast<T>(cellIt.getCorner(p1).getValue());
204
205 // calculate the surface-grid intersection point
206 if (d0 == -d1) { // includes case where d0=d1=0
207 cc[z] = static_cast<T>(cellIt.getIndices(z)) + 0.5;
208 } else {
209 if (std::abs(d0) <= std::abs(d1)) {
210 cc[z] = static_cast<T>(cellIt.getIndices(z)) +
211 (d0 / (d0 - d1));
212 } else {
213 cc[z] = static_cast<T>(cellIt.getIndices(z) + 1) -
214 (d1 / (d1 - d0));
215 }
216 }
217 cc[z] = std::max(
218 cc[z], static_cast<T>(cellIt.getIndices(z) + epsilon));
219 cc[z] = std::min(
220 cc[z],
221 static_cast<T>((cellIt.getIndices(z) + 1) - epsilon));
222 }
223 cc[z] *= gridDelta;
224 }
225
226 int nodeIdx = -1;
227 if (checkNodeFlag) {
228 auto q = quantize(cc);
229 auto it = nodeIdByBin.find(q);
230 if (it != nodeIdByBin.end())
231 nodeIdx = it->second;
232 }
233 if (nodeIdx >= 0) {
234 nod_numbers[n] = nodeIdx;
235 } else {
236 // insert new node
237 nod_numbers[n] =
238 mesh->insertNextNode(cc); // insert new surface node
239 nodes[dir][d] = nod_numbers[n];
240 if (checkNodeFlag)
241 nodeIdByBin.emplace(quantize(cc), nod_numbers[n]);
242
243 for (int a = 0; a < D; a++) {
244 if (cc[a] < mesh->minimumExtent[a])
245 mesh->minimumExtent[a] = cc[a];
246 if (cc[a] > mesh->maximumExtent[a])
247 mesh->maximumExtent[a] = cc[a];
248 }
249 }
250 }
251 }
252
253 // check for degenerate triangles
254 if (!triangleMisformed(nod_numbers)) {
255
256 // only insert if not already present
257 if (uniqueElements.insert(elementI3(nod_numbers)).second) {
258 Vec3D<T> normal;
259 if constexpr (D == 2) {
260 normal = Vec3D<T>{-(mesh->nodes[nod_numbers[1]][1] -
261 mesh->nodes[nod_numbers[0]][1]),
262 mesh->nodes[nod_numbers[1]][0] -
263 mesh->nodes[nod_numbers[0]][0],
264 T(0)};
265 } else {
266 normal = calculateNormal(mesh->nodes[nod_numbers[0]],
267 mesh->nodes[nod_numbers[1]],
268 mesh->nodes[nod_numbers[2]]);
269 }
270
271 double n2 = normal[0] * normal[0] + normal[1] * normal[1] +
272 normal[2] * normal[2];
273 if (n2 > epsilon) {
274 // insert new element
275 mesh->insertNextElement(nod_numbers);
276
277 // insert normal
278 T invn = static_cast<T>(1.) / std::sqrt(static_cast<T>(n2));
279 for (int d = 0; d < D; d++) {
280 normal[d] *= invn;
281 }
282 normals.push_back(normal);
283
284 // insert material
285 if (useMaterialMap) {
286 materials.push_back(materialMap->getMaterialId(l));
287 } else {
288 materials.push_back(static_cast<T>(l));
289 }
290 }
291 }
292 }
293 }
294 }
295 }
296
297 mesh->cellData.insertNextScalarData(materials, "MaterialIds");
298 mesh->cellData.insertNextVectorData(normals, "Normals");
299 mesh->triangles.shrink_to_fit();
300 mesh->nodes.shrink_to_fit();
301 }
302
303private:
304 static inline bool
305 triangleMisformed(const std::array<unsigned, D> &nod_numbers) noexcept {
306 if constexpr (D == 3) {
307 return nod_numbers[0] == nod_numbers[1] ||
308 nod_numbers[0] == nod_numbers[2] ||
309 nod_numbers[1] == nod_numbers[2];
310 } else {
311 return nod_numbers[0] == nod_numbers[1];
312 }
313 }
314
315 static inline Vec3D<NumericType>
316 calculateNormal(const Vec3D<NumericType> &nodeA,
317 const Vec3D<NumericType> &nodeB,
318 const Vec3D<NumericType> &nodeC) noexcept {
319 return CrossProduct(nodeB - nodeA, nodeC - nodeA);
320 }
321};
322
324
325} // 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:79
void setMesh(SmartPointer< viennals::Mesh< T > > passedMesh)
Definition lsToMultiSurfaceMesh.hpp:71
void apply()
Definition lsToMultiSurfaceMesh.hpp:83
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