PolyFEM
Loading...
Searching...
No Matches
RefElementSampler.cpp
Go to the documentation of this file.
3
4#include <igl/predicates/ear_clipping.h>
5
6#include <igl/edges.h>
7#include <igl/boundary_loop.h>
8#include <igl/doublearea.h>
9#include <igl/upsample.h>
10
11#ifdef POLYFEM_WITH_TRIANGLE
12#include <igl/triangle/triangulate.h>
13#endif
14
15#include <cassert>
16
17namespace polyfem
18{
19 using namespace mesh;
20 namespace utils
21 {
22
31 void regular_2d_grid(const int n, bool tri, Eigen::MatrixXd &V, Eigen::MatrixXi &F)
32 {
33
34 V.resize(n * n, 2);
35 F.resize((n - 1) * (n - 1) * 2, 3);
36 const double delta = 1. / (n - 1.);
37 std::vector<int> map(n * n, -1);
38
39 int index = 0;
40 for (int i = 0; i < n; ++i)
41 {
42 for (int j = 0; j < n; ++j)
43 {
44 if (tri && i + j >= n)
45 continue;
46 map[i + j * n] = index;
47 V.row(index) << i * delta, j * delta;
48 ++index;
49 }
50 }
51
52 V.conservativeResize(index, 2);
53
54 std::array<int, 3> tmp;
55
56 index = 0;
57 for (int i = 0; i < n - 1; ++i)
58 {
59 for (int j = 0; j < n - 1; ++j)
60 {
61 tmp = {{map[i + j * n], map[i + 1 + j * n], map[i + (j + 1) * n]}};
62 if (tmp[0] >= 0 && tmp[1] >= 0 && tmp[2] >= 0)
63 {
64 F.row(index) << tmp[0], tmp[1], tmp[2];
65 ++index;
66 }
67
68 tmp = {{map[i + 1 + j * n], map[i + 1 + (j + 1) * n], map[i + (j + 1) * n]}};
69 if (tmp[0] >= 0 && tmp[1] >= 0 && tmp[2] >= 0)
70 {
71 F.row(index) << tmp[0], tmp[1], tmp[2];
72 ++index;
73 }
74 }
75 }
76
77 F.conservativeResize(index, 3);
78 }
79
80 namespace
81 {
82
83 void add_tet(const std::array<int, 4> &tmp, const Eigen::MatrixXd &V, int &index, Eigen::MatrixXi &T)
84 {
85 if (tmp[0] >= 0 && tmp[1] >= 0 && tmp[2] >= 0 && tmp[3] >= 0)
86 {
87 const Eigen::Vector3d e0 = V.row(tmp[1]) - V.row(tmp[0]);
88 const Eigen::Vector3d e1 = V.row(tmp[2]) - V.row(tmp[0]);
89 const Eigen::Vector3d e2 = V.row(tmp[3]) - V.row(tmp[0]);
90 double vol = (e0.cross(e1)).dot(e2);
91 if (vol < 0)
92 T.row(index) << tmp[0], tmp[1], tmp[2], tmp[3];
93 else
94 T.row(index) << tmp[0], tmp[1], tmp[3], tmp[2];
95#ifndef NDEBUG
96 const Eigen::Vector3d ed0 = V.row(T(index, 1)) - V.row(T(index, 0));
97 const Eigen::Vector3d ed1 = V.row(T(index, 2)) - V.row(T(index, 0));
98 const Eigen::Vector3d ed2 = V.row(T(index, 3)) - V.row(T(index, 0));
99 assert((ed0.cross(ed1)).dot(ed2) < 0);
100#endif
101 ++index;
102 }
103 }
104
105 } // anonymous namespace
106
117 void regular_3d_grid(const int nn,
118 bool tet,
119 bool prism,
120 Eigen::MatrixXd &V,
121 Eigen::MatrixXi &F,
122 Eigen::MatrixXi &T)
123 {
124 const int n = nn;
125 const double delta = 1. / (n - 1.);
126 T.resize((n - 1) * (n - 1) * (n - 1) * 6, 4);
127 V.resize(n * n * n, 3);
128 std::vector<int> map(n * n * n, -1);
129
130 int index = 0;
131 for (int i = 0; i < n; ++i)
132 {
133 for (int j = 0; j < n; ++j)
134 {
135 for (int k = 0; k < n; ++k)
136 {
137 if (tet && i + j + k >= n)
138 continue;
139 if (prism && i + j >= n)
140 continue;
141 map[(i + j * n) * n + k] = index;
142 V.row(index) << i * delta, j * delta, k * delta;
143 ++index;
144 }
145 }
146 }
147 V.conservativeResize(index, 3);
148
149 std::array<int, 8> indices;
150 std::array<int, 4> tmp;
151 index = 0;
152 for (int i = 0; i < n - 1; ++i)
153 {
154 for (int j = 0; j < n - 1; ++j)
155 {
156 for (int k = 0; k < n - 1; ++k)
157 {
158 indices = {{(i + j * n) * n + k,
159 (i + 1 + j * n) * n + k,
160 (i + 1 + (j + 1) * n) * n + k,
161 (i + (j + 1) * n) * n + k,
162
163 (i + j * n) * n + k + 1,
164 (i + 1 + j * n) * n + k + 1,
165 (i + 1 + (j + 1) * n) * n + k + 1,
166 (i + (j + 1) * n) * n + k + 1}};
167
168 tmp = {{map[indices[1 - 1]], map[indices[2 - 1]], map[indices[4 - 1]], map[indices[5 - 1]]}};
169 add_tet(tmp, V, index, T);
170
171 tmp = {{map[indices[6 - 1]], map[indices[3 - 1]], map[indices[7 - 1]], map[indices[8 - 1]]}};
172 add_tet(tmp, V, index, T);
173
174 tmp = {{map[indices[5 - 1]], map[indices[2 - 1]], map[indices[6 - 1]], map[indices[4 - 1]]}};
175 add_tet(tmp, V, index, T);
176
177 tmp = {{map[indices[5 - 1]], map[indices[4 - 1]], map[indices[8 - 1]], map[indices[6 - 1]]}};
178 add_tet(tmp, V, index, T);
179
180 tmp = {{map[indices[4 - 1]], map[indices[2 - 1]], map[indices[6 - 1]], map[indices[3 - 1]]}};
181 add_tet(tmp, V, index, T);
182
183 tmp = {{map[indices[3 - 1]], map[indices[4 - 1]], map[indices[8 - 1]], map[indices[6 - 1]]}};
184 add_tet(tmp, V, index, T);
185 }
186 }
187 }
188
189 T.conservativeResize(index, 4);
190
191 F.resize(4 * index, 3);
192
193 F.block(0, 0, index, 1) = T.col(1);
194 F.block(0, 1, index, 1) = T.col(0);
195 F.block(0, 2, index, 1) = T.col(2);
196
197 F.block(index, 0, index, 1) = T.col(0);
198 F.block(index, 1, index, 1) = T.col(1);
199 F.block(index, 2, index, 1) = T.col(3);
200
201 F.block(2 * index, 0, index, 1) = T.col(1);
202 F.block(2 * index, 1, index, 1) = T.col(2);
203 F.block(2 * index, 2, index, 1) = T.col(3);
204
205 F.block(3 * index, 0, index, 1) = T.col(2);
206 F.block(3 * index, 1, index, 1) = T.col(0);
207 F.block(3 * index, 2, index, 1) = T.col(3);
208 }
209
210 void RefElementSampler::init(const bool is_volume, const int n_elements, double target_rel_area)
211 {
212 is_volume_ = is_volume;
213
214 area_param_ = target_rel_area * n_elements;
215#ifndef NDEBUG
216 area_param_ *= 10.0;
217#endif
218
219 build();
220 }
221
223 {
224 using namespace Eigen;
225
226 if (is_volume_)
227 {
228 // cube
229 {
230 MatrixXd pts(8, 3);
231 pts << 0, 0, 0,
232 0, 1, 0,
233 1, 1, 0,
234 1, 0, 0,
235
236 // 4
237 0, 0, 1,
238 0, 1, 1,
239 1, 1, 1,
240 1, 0, 1;
241
242 Eigen::MatrixXi faces(12, 3);
243 faces << 1, 2, 0,
244 0, 2, 3,
245
246 5, 4, 6,
247 4, 7, 6,
248
249 1, 0, 4,
250 1, 4, 5,
251
252 2, 1, 5,
253 2, 5, 6,
254
255 3, 2, 6,
256 3, 6, 7,
257
258 0, 3, 7,
259 0, 7, 4;
260
261 regular_3d_grid(std::max(2., round(1. / pow(area_param_, 1. / 3.) + 1) / 2.), false, false, cube_points_, cube_faces_, cube_tets_);
262
263 // Extract sampled edges matching the base element edges
264 Eigen::MatrixXi edges(12, 2);
265 edges << 0, 1,
266 1, 2,
267 2, 3,
268 3, 0,
269
270 4, 5,
271 5, 6,
272 6, 7,
273 7, 4,
274
275 0, 4,
276 1, 5,
277 2, 6,
278 3, 7;
279 igl::edges(cube_faces_, cube_edges_);
281
282 // Same local order as in FEMBasis3d
283 cube_corners_.resize(8, 3);
284 cube_corners_ << 0, 0, 0,
285 1, 0, 0,
286 1, 1, 0,
287 0, 1, 0,
288 0, 0, 1,
289 1, 0, 1,
290 1, 1, 1,
291 0, 1, 1;
292 }
293
294 // tet
295 {
296 MatrixXd pts(4, 3);
297 pts << 0, 0, 0,
298 1, 0, 0,
299 0, 1, 0,
300 0, 0, 1;
301
302 Eigen::MatrixXi faces(4, 3);
303 faces << 0, 1, 2,
304
305 3, 1, 0,
306 2, 1, 3,
307 0, 2, 3;
308
309 regular_3d_grid(std::max(2., round(1. / pow(area_param_, 1. / 3.) + 1)), true, false, simplex_points_, simplex_faces_, simplex_tets_);
310
311 // Extract sampled edges matching the base element edges
312 Eigen::MatrixXi edges;
313 igl::edges(faces, edges);
314 igl::edges(simplex_faces_, simplex_edges_);
316
317 // Same local order as in FEMBasis3d
318 simplex_corners_.resize(4, 3);
319 simplex_corners_ << 0, 0, 0,
320 1, 0, 0,
321 0, 1, 0,
322 0, 0, 1;
323 }
324
325 // prism
326 {
327 MatrixXd pts(6, 3);
328 pts << 0, 0, 0,
329 0, 1, 0,
330 1, 0, 0,
331 0, 0, 1,
332 0, 1, 1,
333 1, 0, 1;
334
335 regular_3d_grid(std::max(2., round(1. / pow(area_param_, 1. / 3.) + 1) / 2.),
336 false, true, prism_points_, prism_faces_, prism_tets_);
337
338 // Extract sampled edges matching the base element edges
339 Eigen::MatrixXi edges(9, 2);
340 edges << 0, 1,
341 1, 2,
342 2, 0,
343
344 3, 4,
345 4, 5,
346 5, 3,
347
348 0, 3,
349 1, 4,
350 2, 5;
351
352 igl::edges(prism_faces_, prism_edges_);
354
355 // Same local order as in FEMBasis3d
356 prism_corners_.resize(6, 3);
357 prism_corners_ << 0, 0, 0,
358 0, 1, 0,
359 1, 0, 0,
360 0, 0, 1,
361 0, 1, 1,
362 1, 0, 1;
363 }
364 }
365 else
366 {
367 {
368 MatrixXd pts(4, 2);
369 pts << 0, 0,
370 0, 1,
371 1, 1,
372 1, 0;
373
374 MatrixXi E(4, 2);
375 E << 0, 1,
376 1, 2,
377 2, 3,
378 3, 0;
379
380 regular_2d_grid(std::max(2., round(1. / sqrt(area_param_) + 1)), false, cube_points_, cube_faces_);
381
382 // Extract sampled edges matching the base element edges
383 igl::edges(cube_faces_, cube_edges_);
385
386 // Same local order as in FEMBasis2d
387 cube_corners_.resize(4, 2);
388 cube_corners_ << 0, 0,
389 1, 0,
390 1, 1,
391 0, 1;
392 }
393 {
394 MatrixXd pts(3, 2);
395 pts << 0, 0,
396 0, 1,
397 1, 0;
398
399 MatrixXi E(3, 2);
400 E << 0, 1,
401 1, 2,
402 2, 0;
403
404 regular_2d_grid(std::max(2., round(1. / sqrt(area_param_) + 1)), true, simplex_points_, simplex_faces_);
405
406 // Extract sampled edges matching the base element edges
407 igl::edges(simplex_faces_, simplex_edges_);
409
410 // Same local order as in FEMBasis2d
411 simplex_corners_.resize(3, 2);
412 simplex_corners_ << 0, 0,
413 1, 0,
414 0, 1;
415 }
416 }
417 }
418
419 void RefElementSampler::sample_polygon(const Eigen::MatrixXd &poly, Eigen::MatrixXd &pts, Eigen::MatrixXi &faces, Eigen::MatrixXi &edges) const
420 {
421
422#ifdef POLYFEM_WITH_TRIANGLE
423 Eigen::MatrixXi E(poly.rows(), 2);
424 const Eigen::MatrixXd H(0, 2);
425 const std::string flags = "Qzqa" + std::to_string(area_param_ / 10.0);
426
427 for (int i = 0; i < poly.rows(); ++i)
428 E.row(i) << i, (i + 1) % poly.rows();
429
430 igl::triangle::triangulate(poly, E, H, flags, pts, faces);
431#else
432
433 const Eigen::MatrixXi rt = Eigen::MatrixXi::Zero(poly.rows(), 1);
434 faces.resize(0, 0);
435 Eigen::VectorXi I;
436 Eigen::MatrixXd area;
437 igl::predicates::ear_clipping(poly, rt, faces, I);
438
439 igl::doublearea(poly, faces, area);
440
441 const double area_avg = area.array().sum() / poly.rows() / 2;
442
443 const int n_refs = area_avg / area_param_ * 40;
444
445 Eigen::MatrixXi new_faces;
446
447 igl::upsample(poly, faces, pts, new_faces, n_refs);
448
449 faces = new_faces;
450#endif
451 std::vector<int> loop;
452 igl::default_num_threads(1);
453 igl::boundary_loop(faces, loop);
454 igl::default_num_threads(0);
455 edges.resize(loop.size(), 2);
456 for (int i = 0; i < loop.size(); ++i)
457 edges.row(i) << loop[i], loop[(i + 1) % loop.size()];
458 }
459
460 void RefElementSampler::sample_polyhedron(const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &f, Eigen::MatrixXd &pts, Eigen::MatrixXi &tets, Eigen::MatrixXi &faces) const
461 {
462 const Eigen::MatrixXd kernel = vertices.colwise().mean();
463
464 polyfem::mesh::tertrahedralize_star_shaped_surface(vertices, f, kernel, pts, faces, tets);
465 }
466 } // namespace utils
467} // namespace polyfem
int V
std::vector< Eigen::VectorXi > faces
void init(const bool is_volume, const int n_elements, const double target_rel_area)
void sample_polygon(const Eigen::MatrixXd &poly, Eigen::MatrixXd &pts, Eigen::MatrixXi &faces, Eigen::MatrixXi &edges) const
void sample_polyhedron(const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &f, Eigen::MatrixXd &pts, Eigen::MatrixXi &faces, Eigen::MatrixXi &edges) const
void extract_parent_edges(const Eigen::MatrixXd &IV, const Eigen::MatrixXi &IE, const Eigen::MatrixXd &BV, const Eigen::MatrixXi &BE, Eigen::MatrixXi &OE)
Extract a set of edges that are overlap with a set given set of parent edges, using vertices position...
void tertrahedralize_star_shaped_surface(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::RowVector3d &kernel, Eigen::MatrixXd &OV, Eigen::MatrixXi &OF, Eigen::MatrixXi &OT)
Tetrahedralize a star-shaped mesh, with a given point in its kernel.
void regular_2d_grid(const int n, bool tri, Eigen::MatrixXd &V, Eigen::MatrixXi &F)
Generate a canonical triangle/quad subdivided from a regular grid.
void regular_3d_grid(const int nn, bool tet, bool prism, Eigen::MatrixXd &V, Eigen::MatrixXi &F, Eigen::MatrixXi &T)
Generate a canonical tet/hex subdivided from a regular grid.