PolyFEM
Loading...
Searching...
No Matches
BoundarySampler.cpp
Go to the documentation of this file.
1#include "BoundarySampler.hpp"
2
6
10
13
14#include <cassert>
15
16namespace polyfem
17{
18 using namespace mesh;
19 using namespace quadrature;
20 namespace utils
21 {
22 namespace
23 {
24 Eigen::RowVector2d quad_local_node_coordinates(int local_index)
25 {
26 assert(local_index >= 0 && local_index < 4);
27 Eigen::MatrixXd p;
29 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
30 }
31
32 Eigen::RowVector2d tri_local_node_coordinates(int local_index)
33 {
34 Eigen::MatrixXd p;
36 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
37 }
38
39 Eigen::RowVector3d linear_tet_local_node_coordinates(int local_index)
40 {
41 Eigen::MatrixXd p;
43 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
44 }
45
46 Eigen::RowVector3d linear_hex_local_node_coordinates(int local_index)
47 {
48 Eigen::MatrixXd p;
50 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
51 }
52
53 Eigen::RowVector3d linear_prism_local_node_coordinates(int local_index)
54 {
55 Eigen::MatrixXd p;
57 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
58 }
59
60 Eigen::Matrix2d quad_local_node_coordinates_from_edge(int le)
61 {
62 Eigen::Matrix2d res(2, 2);
63 res.row(0) = quad_local_node_coordinates(le);
64 res.row(1) = quad_local_node_coordinates((le + 1) % 4);
65
66 return res;
67 }
68
69 Eigen::Matrix2d tri_local_node_coordinates_from_edge(int le)
70 {
71 Eigen::Matrix2d res(2, 2);
72 res.row(0) = tri_local_node_coordinates(le);
73 res.row(1) = tri_local_node_coordinates((le + 1) % 3);
74
75 return res;
76 }
77 } // namespace
78
80 {
81 Eigen::Matrix2d res(2, 2);
82 res.row(0) = quad_local_node_coordinates(le);
83 res.row(1) = quad_local_node_coordinates((le + 1) % 4);
84
85 return res;
86 }
87
89 {
90 Eigen::Matrix2d res(2, 2);
91 res.row(0) = tri_local_node_coordinates(le);
92 res.row(1) = tri_local_node_coordinates((le + 1) % 3);
93
94 return res;
95 }
96
98 {
99 Eigen::Matrix<int, 4, 3> fv;
100 fv.row(0) << 0, 1, 2;
101 fv.row(1) << 0, 1, 3;
102 fv.row(2) << 1, 2, 3;
103 fv.row(3) << 2, 0, 3;
104
105 Eigen::MatrixXd res(3, 3);
106 for (int i = 0; i < 3; ++i)
107 res.row(i) = linear_tet_local_node_coordinates(fv(lf, i));
108
109 return res;
110 }
111
113 {
114 Eigen::Matrix<int, 6, 4> fv;
115 fv.row(0) << 0, 3, 7, 4;
116 fv.row(1) << 1, 2, 6, 5;
117 fv.row(2) << 0, 1, 5, 4;
118 fv.row(3) << 3, 2, 6, 7;
119 fv.row(4) << 0, 1, 2, 3;
120 fv.row(5) << 4, 5, 6, 7;
121
122 Eigen::MatrixXd res(4, 3);
123 for (int i = 0; i < 4; ++i)
124 res.row(i) = linear_hex_local_node_coordinates(fv(lf, i));
125
126 return res;
127 }
128
130 {
131 Eigen::Matrix<int, 5, 4> fv;
132 fv.row(0) << 0, 1, 2, -1;
133 fv.row(1) << 3, 4, 5, -1;
134
135 fv.row(2) << 0, 1, 4, 3;
136 fv.row(3) << 1, 2, 5, 4;
137 fv.row(4) << 2, 0, 3, 5;
138
139 const int nv = lf < 2 ? 3 : 4;
140 Eigen::MatrixXd res(nv, 3);
141 for (int i = 0; i < nv; ++i)
142 res.row(i) = linear_prism_local_node_coordinates(fv(lf, i));
143
144 return res;
145 }
146
147 void utils::BoundarySampler::quadrature_for_quad_edge(int index, int order, int gid, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
148 {
149 auto endpoints = quad_local_node_coordinates_from_edge(index);
150
151 Quadrature quad;
152 LineQuadrature quad_rule;
153 quad_rule.get_quadrature(order, quad);
154
155 points.resize(quad.points.rows(), endpoints.cols());
156 uv.resize(quad.points.rows(), 2);
157
158 uv.col(0) = (1.0 - quad.points.array());
159 uv.col(1) = quad.points.array();
160
161 for (int c = 0; c < 2; ++c)
162 {
163 points.col(c) = (1.0 - quad.points.array()) * endpoints(0, c) + quad.points.array() * endpoints(1, c);
164 }
165
166 weights = quad.weights;
167 weights *= mesh.edge_length(gid);
168 }
169
170 void utils::BoundarySampler::quadrature_for_tri_edge(int index, int order, int gid, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
171 {
172 auto endpoints = tri_local_node_coordinates_from_edge(index);
173
174 Quadrature quad;
175 LineQuadrature quad_rule;
176 quad_rule.get_quadrature(order, quad);
177
178 points.resize(quad.points.rows(), endpoints.cols());
179 uv.resize(quad.points.rows(), 2);
180
181 uv.col(0) = (1.0 - quad.points.array());
182 uv.col(1) = quad.points.array();
183
184 for (int c = 0; c < 2; ++c)
185 {
186 points.col(c) = (1.0 - quad.points.array()) * endpoints(0, c) + quad.points.array() * endpoints(1, c);
187 }
188
189 weights = quad.weights;
190 weights *= mesh.edge_length(gid);
191 }
192
193 void utils::BoundarySampler::quadrature_for_quad_face(int index, int order, int gid, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
194 {
195 auto endpoints = hex_local_node_coordinates_from_face(index);
196
197 Quadrature quad;
198 QuadQuadrature quad_rule;
199 quad_rule.get_quadrature(order, quad);
200
201 const int n_pts = quad.points.rows();
202 points.resize(n_pts, endpoints.cols());
203
204 uv.resize(quad.points.rows(), 4);
205 uv.col(0) = quad.points.col(0);
206 uv.col(1) = 1 - uv.col(0).array();
207 uv.col(2) = quad.points.col(1);
208 uv.col(3) = 1 - uv.col(2).array();
209
210 for (int i = 0; i < n_pts; ++i)
211 {
212 const double b1 = quad.points(i, 0);
213 const double b2 = 1 - b1;
214
215 const double b3 = quad.points(i, 1);
216 const double b4 = 1 - b3;
217
218 for (int c = 0; c < 3; ++c)
219 {
220 points(i, c) = b3 * (b1 * endpoints(0, c) + b2 * endpoints(1, c)) + b4 * (b1 * endpoints(3, c) + b2 * endpoints(2, c));
221 }
222 }
223
224 weights = quad.weights;
225 weights *= mesh.quad_area(gid);
226 }
227
228 void utils::BoundarySampler::quadrature_for_tri_face(int index, int order, int gid, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
229 {
230 auto endpoints = tet_local_node_coordinates_from_face(index);
231 Quadrature quad;
232 TriQuadrature quad_rule;
233 quad_rule.get_quadrature(order, quad);
234
235 const int n_pts = quad.points.rows();
236 points.resize(n_pts, endpoints.cols());
237
238 uv.resize(quad.points.rows(), 3);
239 uv.col(0) = quad.points.col(0);
240 uv.col(1) = quad.points.col(1);
241 uv.col(2) = 1 - uv.col(0).array() - uv.col(1).array();
242
243 for (int i = 0; i < n_pts; ++i)
244 {
245 const double b1 = quad.points(i, 0);
246 const double b3 = quad.points(i, 1);
247 const double b2 = 1 - b1 - b3;
248
249 for (int c = 0; c < 3; ++c)
250 {
251 points(i, c) = b1 * endpoints(0, c) + b2 * endpoints(1, c) + b3 * endpoints(2, c);
252 }
253 }
254
255 weights = quad.weights;
256 // 2 * because weights sum to 1/2 already
257 weights *= 2 * mesh.tri_area(gid);
258 }
259
260 void utils::BoundarySampler::quadrature_for_prism_face(int index, int orderp, int orderq, int gid, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
261 {
262 auto endpoints = prism_local_node_coordinates_from_face(index);
263
264 if (index < 2)
265 {
266 Quadrature quad;
267 TriQuadrature quad_rule;
268 quad_rule.get_quadrature(orderp, quad);
269
270 const int n_pts = quad.points.rows();
271 points.resize(n_pts, endpoints.cols());
272
273 uv.resize(quad.points.rows(), 3);
274 uv.col(0) = quad.points.col(0);
275 uv.col(1) = quad.points.col(1);
276 uv.col(2) = 1 - uv.col(0).array() - uv.col(1).array();
277
278 for (int i = 0; i < n_pts; ++i)
279 {
280 const double b1 = quad.points(i, 0);
281 const double b3 = quad.points(i, 1);
282 const double b2 = 1 - b1 - b3;
283
284 for (int c = 0; c < 3; ++c)
285 {
286 points(i, c) = b1 * endpoints(0, c) + b2 * endpoints(1, c) + b3 * endpoints(2, c);
287 }
288 }
289
290 weights = quad.weights;
291 // 2 * because weights sum to 1/2 already
292 weights *= 2 * mesh.tri_area(gid);
293 }
294 else
295 {
296 Quadrature quad;
297 QuadQuadrature quad_rule;
298 quad_rule.get_quadrature(orderq, quad);
299
300 const int n_pts = quad.points.rows();
301 points.resize(n_pts, endpoints.cols());
302
303 uv.resize(quad.points.rows(), 4);
304 uv.col(0) = quad.points.col(0);
305 uv.col(1) = 1 - uv.col(0).array();
306 uv.col(2) = quad.points.col(1);
307 uv.col(3) = 1 - uv.col(2).array();
308
309 for (int i = 0; i < n_pts; ++i)
310 {
311 const double b1 = quad.points(i, 0);
312 const double b2 = 1 - b1;
313
314 const double b3 = quad.points(i, 1);
315 const double b4 = 1 - b3;
316
317 for (int c = 0; c < 3; ++c)
318 {
319 points(i, c) = b3 * (b1 * endpoints(0, c) + b2 * endpoints(1, c)) + b4 * (b1 * endpoints(3, c) + b2 * endpoints(2, c));
320 }
321 }
322
323 weights = quad.weights;
324 weights *= mesh.quad_area(gid);
325 }
326 }
327
328 void utils::BoundarySampler::sample_parametric_quad_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
329 {
330 auto endpoints = quad_local_node_coordinates_from_edge(index);
331 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
332 samples.resize(n_samples, endpoints.cols());
333
334 uv.resize(n_samples, 2);
335
336 uv.col(0) = (1.0 - t.array());
337 uv.col(1) = t.array();
338
339 for (int c = 0; c < 2; ++c)
340 {
341 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
342 }
343 }
344
345 void utils::BoundarySampler::sample_parametric_tri_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
346 {
347 auto endpoints = tri_local_node_coordinates_from_edge(index);
348 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
349 samples.resize(n_samples, endpoints.cols());
350
351 uv.resize(n_samples, 2);
352
353 uv.col(0) = (1.0 - t.array());
354 uv.col(1) = t.array();
355
356 for (int c = 0; c < 2; ++c)
357 {
358 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
359 }
360 }
361
362 void utils::BoundarySampler::sample_parametric_quad_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
363 {
364 auto endpoints = hex_local_node_coordinates_from_face(index);
365 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
366 samples.resize(n_samples * n_samples, endpoints.cols());
367 Eigen::MatrixXd left(n_samples, endpoints.cols());
368 Eigen::MatrixXd right(n_samples, endpoints.cols());
369
370 uv.resize(n_samples * n_samples, endpoints.cols());
371 uv.setZero();
372
373 for (int c = 0; c < 3; ++c)
374 {
375 left.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
376 right.col(c) = (1.0 - t.array()).matrix() * endpoints(3, c) + t * endpoints(2, c);
377 }
378 for (int c = 0; c < 3; ++c)
379 {
380 Eigen::MatrixXd x = (1.0 - t.array()).matrix() * left.col(c).transpose() + t * right.col(c).transpose();
381 assert(x.size() == n_samples * n_samples);
382
383 samples.col(c) = Eigen::Map<const Eigen::VectorXd>(x.data(), x.size());
384 }
385 }
386
387 void utils::BoundarySampler::sample_parametric_tri_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
388 {
389 auto endpoints = tet_local_node_coordinates_from_face(index);
390 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
391 samples.resize(n_samples * n_samples, endpoints.cols());
392
393 uv.resize(n_samples * n_samples, 3);
394
395 int counter = 0;
396 for (int u = 0; u < n_samples; ++u)
397 {
398 for (int v = 0; v < n_samples; ++v)
399 {
400 if (t(u) + t(v) > 1)
401 continue;
402
403 uv(counter, 0) = t(u);
404 uv(counter, 1) = t(v);
405 uv(counter, 2) = 1 - t(u) - t(v);
406
407 for (int c = 0; c < 3; ++c)
408 {
409 samples(counter, c) = t(u) * endpoints(0, c) + t(v) * endpoints(1, c) + (1 - t(u) - t(v)) * endpoints(2, c);
410 }
411 ++counter;
412 }
413 }
414 samples.conservativeResize(counter, 3);
415 uv.conservativeResize(counter, 3);
416 }
417
418 void utils::BoundarySampler::sample_parametric_prism_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
419 {
420 auto endpoints = prism_local_node_coordinates_from_face(index);
421
422 if (index < 2)
423 {
424 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
425 samples.resize(n_samples * n_samples, endpoints.cols());
426
427 uv.resize(n_samples * n_samples, 3);
428
429 int counter = 0;
430 for (int u = 0; u < n_samples; ++u)
431 {
432 for (int v = 0; v < n_samples; ++v)
433 {
434 if (t(u) + t(v) > 1)
435 continue;
436
437 uv(counter, 0) = t(u);
438 uv(counter, 1) = t(v);
439 uv(counter, 2) = 1 - t(u) - t(v);
440
441 for (int c = 0; c < 3; ++c)
442 {
443 samples(counter, c) = t(u) * endpoints(0, c) + t(v) * endpoints(1, c) + (1 - t(u) - t(v)) * endpoints(2, c);
444 }
445 ++counter;
446 }
447 }
448 samples.conservativeResize(counter, 3);
449 uv.conservativeResize(counter, 3);
450 }
451 else
452 {
453 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
454 samples.resize(n_samples * n_samples, endpoints.cols());
455 Eigen::MatrixXd left(n_samples, endpoints.cols());
456 Eigen::MatrixXd right(n_samples, endpoints.cols());
457
458 uv.resize(n_samples * n_samples, endpoints.cols());
459 uv.setZero();
460
461 for (int c = 0; c < 3; ++c)
462 {
463 left.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
464 right.col(c) = (1.0 - t.array()).matrix() * endpoints(3, c) + t * endpoints(2, c);
465 }
466 for (int c = 0; c < 3; ++c)
467 {
468 Eigen::MatrixXd x = (1.0 - t.array()).matrix() * left.col(c).transpose() + t * right.col(c).transpose();
469 assert(x.size() == n_samples * n_samples);
470
471 samples.col(c) = Eigen::Map<const Eigen::VectorXd>(x.data(), x.size());
472 }
473 }
474 }
475
476 void utils::BoundarySampler::sample_polygon_edge(int face_id, int edge_id, int n_samples, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
477 {
478 assert(!mesh.is_volume());
479
480 const CMesh2D &mesh2d = dynamic_cast<const CMesh2D &>(mesh);
481
482 auto index = mesh2d.get_index_from_face(face_id);
483
484 bool found = false;
485 for (int i = 0; i < mesh2d.n_face_vertices(face_id); ++i)
486 {
487 if (index.edge == edge_id)
488 {
489 found = true;
490 break;
491 }
492
493 index = mesh2d.next_around_face(index);
494 }
495
496 assert(found);
497
498 auto p0 = mesh2d.point(index.vertex);
499 auto p1 = mesh2d.point(mesh2d.switch_edge(index).vertex);
500 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
501 samples.resize(n_samples, p0.cols());
502
503 uv.resize(n_samples, 2);
504
505 uv.col(0) = (1.0 - t.array());
506 uv.col(1) = t.array();
507
508 for (int c = 0; c < 2; ++c)
509 {
510 samples.col(c) = (1.0 - t.array()) * p0(c) + t.array() * p1(c);
511 }
512 }
513
514 void utils::BoundarySampler::quadrature_for_polygon_edge(int face_id, int edge_id, int order, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
515 {
516 assert(!mesh.is_volume());
517
518 const CMesh2D &mesh2d = dynamic_cast<const CMesh2D &>(mesh);
519
520 auto index = mesh2d.get_index_from_face(face_id);
521
522 bool found = false;
523 for (int i = 0; i < mesh2d.n_face_vertices(face_id); ++i)
524 {
525 if (index.edge == edge_id)
526 {
527 found = true;
528 break;
529 }
530
531 index = mesh2d.next_around_face(index);
532 }
533
534 assert(found);
535
536 auto p0 = mesh2d.point(index.vertex);
537 auto p1 = mesh2d.point(mesh2d.switch_edge(index).vertex);
538 Quadrature quad;
539 LineQuadrature quad_rule;
540 quad_rule.get_quadrature(order, quad);
541
542 points.resize(quad.points.rows(), p0.cols());
543 uv.resize(quad.points.rows(), 2);
544
545 uv.col(0) = (1.0 - quad.points.array());
546 uv.col(1) = quad.points.array();
547
548 for (int c = 0; c < 2; ++c)
549 {
550 points.col(c) = (1.0 - quad.points.array()) * p0(c) + quad.points.array() * p1(c);
551 }
552
553 weights = quad.weights;
554 weights *= mesh.edge_length(edge_id);
555 }
556
557 void utils::BoundarySampler::normal_for_quad_edge(int index, Eigen::MatrixXd &normal)
558 {
559 auto endpoints = quad_local_node_coordinates_from_edge(index);
560 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
561 normal.resize(1, 2);
562 normal(0) = -e(1);
563 normal(1) = e(0);
564 normal.normalize();
565 }
566
567 void utils::BoundarySampler::normal_for_tri_edge(int index, Eigen::MatrixXd &normal)
568 {
569 auto endpoints = tri_local_node_coordinates_from_edge(index);
570 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
571 normal.resize(1, 2);
572 normal(0) = -e(1);
573 normal(1) = e(0);
574 normal.normalize();
575 }
576
577 void utils::BoundarySampler::normal_for_quad_face(int index, Eigen::MatrixXd &normal)
578 {
579 auto endpoints = hex_local_node_coordinates_from_face(index);
580 const Eigen::RowVector3d e1 = endpoints.row(0) - endpoints.row(1);
581 const Eigen::RowVector3d e2 = endpoints.row(0) - endpoints.row(2);
582 normal = e1.cross(e2);
583 normal.normalize();
584
585 if (index == 0 || index == 3 || index == 4)
586 normal *= -1;
587 }
588
589 void utils::BoundarySampler::normal_for_tri_face(int index, Eigen::MatrixXd &normal)
590 {
591 auto endpoints = tet_local_node_coordinates_from_face(index);
592 const Eigen::RowVector3d e1 = endpoints.row(0) - endpoints.row(1);
593 const Eigen::RowVector3d e2 = endpoints.row(0) - endpoints.row(2);
594 normal = e1.cross(e2);
595 normal.normalize();
596
597 if (index == 0)
598 normal *= -1;
599 }
600
601 void utils::BoundarySampler::normal_for_prism_face(int index, Eigen::MatrixXd &normal)
602 {
603 auto endpoints = prism_local_node_coordinates_from_face(index);
604 const Eigen::RowVector3d e1 = endpoints.row(0) - endpoints.row(1);
605 const Eigen::RowVector3d e2 = endpoints.row(0) - endpoints.row(2);
606 normal = e1.cross(e2);
607 normal.normalize();
608
609 if (index == 0)
610 normal *= -1;
611 }
612
613 void utils::BoundarySampler::normal_for_polygon_edge(int face_id, int edge_id, const Mesh &mesh, Eigen::MatrixXd &normal)
614 {
615 assert(!mesh.is_volume());
616
617 const CMesh2D &mesh2d = dynamic_cast<const CMesh2D &>(mesh);
618
619 auto index = mesh2d.get_index_from_face(face_id);
620
621 bool found = false;
622 for (int i = 0; i < mesh2d.n_face_vertices(face_id); ++i)
623 {
624 if (index.edge == edge_id)
625 {
626 found = true;
627 break;
628 }
629
630 index = mesh2d.next_around_face(index);
631 }
632
633 assert(found);
634
635 auto p0 = mesh2d.point(index.vertex);
636 auto p1 = mesh2d.point(mesh2d.switch_edge(index).vertex);
637 const Eigen::MatrixXd e = p0 - p1;
638 normal.resize(1, 2);
639 normal(0) = -e(1);
640 normal(1) = e(0);
641 normal.normalize();
642 }
643
644 bool utils::BoundarySampler::boundary_quadrature(const LocalBoundary &local_boundary, const QuadratureOrders &order, const Mesh &mesh, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::MatrixXd &normals, Eigen::VectorXd &weights, Eigen::VectorXi &global_primitive_ids)
645 {
646 uv.resize(0, 0);
647 points.resize(0, 0);
648 normals.resize(0, 0);
649 weights.resize(0);
650 global_primitive_ids.resize(0);
651
652 for (int i = 0; i < local_boundary.size(); ++i)
653 {
654 const int gid = local_boundary.global_primitive_id(i);
655 Eigen::MatrixXd tmp_p, tmp_uv, tmp_n;
656 Eigen::VectorXd tmp_w;
657 switch (local_boundary.type())
658 {
659 case BoundaryType::TRI_LINE:
660 quadrature_for_tri_edge(local_boundary[i], order[0], gid, mesh, tmp_uv, tmp_p, tmp_w);
661 normal_for_tri_edge(local_boundary[i], tmp_n);
662 break;
663 case BoundaryType::QUAD_LINE:
664 quadrature_for_quad_edge(local_boundary[i], order[0], gid, mesh, tmp_uv, tmp_p, tmp_w);
665 normal_for_quad_edge(local_boundary[i], tmp_n);
666 break;
667 case BoundaryType::QUAD:
668 quadrature_for_quad_face(local_boundary[i], order[0], gid, mesh, tmp_uv, tmp_p, tmp_w);
669 normal_for_quad_face(local_boundary[i], tmp_n);
670 break;
671 case BoundaryType::TRI:
672 quadrature_for_tri_face(local_boundary[i], order[0], gid, mesh, tmp_uv, tmp_p, tmp_w);
673 normal_for_tri_face(local_boundary[i], tmp_n);
674 break;
675 case BoundaryType::PRISM:
676 quadrature_for_prism_face(local_boundary[i], order[0], order[1], gid, mesh, tmp_uv, tmp_p, tmp_w);
677 normal_for_prism_face(local_boundary[i], tmp_n);
678 break;
679 case BoundaryType::POLYGON:
680 quadrature_for_polygon_edge(local_boundary.element_id(), local_boundary.global_primitive_id(i), order[0], mesh, tmp_uv, tmp_p, tmp_w);
681 normal_for_polygon_edge(local_boundary.element_id(), local_boundary.global_primitive_id(i), mesh, tmp_n);
682 break;
683 case BoundaryType::INVALID:
684 assert(false);
685 break;
686 default:
687 assert(false);
688 }
689
690 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
691 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
692
693 points.conservativeResize(points.rows() + tmp_p.rows(), tmp_p.cols());
694 points.bottomRows(tmp_p.rows()) = tmp_p;
695
696 normals.conservativeResize(normals.rows() + tmp_p.rows(), tmp_p.cols());
697 for (int k = normals.rows() - tmp_p.rows(); k < normals.rows(); ++k)
698 normals.row(k) = tmp_n;
699
700 weights.conservativeResize(weights.rows() + tmp_w.rows(), tmp_w.cols());
701 weights.bottomRows(tmp_w.rows()) = tmp_w;
702
703 global_primitive_ids.conservativeResize(global_primitive_ids.rows() + tmp_p.rows());
704 global_primitive_ids.bottomRows(tmp_p.rows()).setConstant(gid);
705 }
706
707 assert(uv.rows() == global_primitive_ids.size());
708 assert(points.rows() == global_primitive_ids.size());
709 assert(normals.rows() == global_primitive_ids.size());
710 assert(weights.size() == global_primitive_ids.size());
711
712 return true;
713 }
714
715 bool utils::BoundarySampler::sample_boundary(const LocalBoundary &local_boundary, const int n_samples, const Mesh &mesh, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples, Eigen::VectorXi &global_primitive_ids)
716 {
717 uv.resize(0, 0);
718 samples.resize(0, 0);
719 global_primitive_ids.resize(0);
720
721 for (int i = 0; i < local_boundary.size(); ++i)
722 {
723 Eigen::MatrixXd tmp, tmp_uv;
724 switch (local_boundary.type())
725 {
726 case BoundaryType::TRI_LINE:
727 sample_parametric_tri_edge(local_boundary[i], n_samples, tmp_uv, tmp);
728 break;
729 case BoundaryType::QUAD_LINE:
730 sample_parametric_quad_edge(local_boundary[i], n_samples, tmp_uv, tmp);
731 break;
732 case BoundaryType::QUAD:
733 sample_parametric_quad_face(local_boundary[i], n_samples, tmp_uv, tmp);
734 break;
735 case BoundaryType::TRI:
736 sample_parametric_tri_face(local_boundary[i], n_samples, tmp_uv, tmp);
737 break;
738 case BoundaryType::PRISM:
739 sample_parametric_prism_face(local_boundary[i], n_samples, tmp_uv, tmp);
740 break;
741 case BoundaryType::POLYGON:
742 sample_polygon_edge(local_boundary.element_id(), local_boundary.global_primitive_id(i), n_samples, mesh, tmp_uv, tmp);
743 break;
744 case BoundaryType::INVALID:
745 assert(false);
746 break;
747 default:
748 assert(false);
749 }
750
751 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
752 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
753
754 samples.conservativeResize(samples.rows() + tmp.rows(), tmp.cols());
755 samples.bottomRows(tmp.rows()) = tmp;
756
757 global_primitive_ids.conservativeResize(global_primitive_ids.rows() + tmp.rows());
758 global_primitive_ids.bottomRows(tmp.rows()).setConstant(local_boundary.global_primitive_id(i));
759 }
760
761 assert(uv.rows() == global_primitive_ids.size());
762 assert(samples.rows() == global_primitive_ids.size());
763
764 return true;
765 }
766
767 bool utils::BoundarySampler::boundary_quadrature(const mesh::LocalBoundary &local_boundary, const QuadratureOrders &order, const mesh::Mesh &mesh, const int i, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::MatrixXd &normals, Eigen::VectorXd &weights)
768 {
769 assert(local_boundary.size() > i);
770
771 uv.resize(0, 0);
772 points.resize(0, 0);
773 weights.resize(0);
774 const int gid = local_boundary.global_primitive_id(i);
775
776 Eigen::MatrixXd normal;
777 switch (local_boundary.type())
778 {
779 case BoundaryType::TRI_LINE:
780 quadrature_for_tri_edge(local_boundary[i], order[0], gid, mesh, uv, points, weights);
781 normal_for_tri_edge(local_boundary[i], normal);
782 break;
783 case BoundaryType::QUAD_LINE:
784 quadrature_for_quad_edge(local_boundary[i], order[0], gid, mesh, uv, points, weights);
785 normal_for_quad_edge(local_boundary[i], normal);
786 break;
787 case BoundaryType::QUAD:
788 quadrature_for_quad_face(local_boundary[i], order[0], gid, mesh, uv, points, weights);
789 normal_for_quad_face(local_boundary[i], normal);
790 break;
791 case BoundaryType::TRI:
792 quadrature_for_tri_face(local_boundary[i], order[0], gid, mesh, uv, points, weights);
793 normal_for_tri_face(local_boundary[i], normal);
794 break;
795 case BoundaryType::PRISM:
796 quadrature_for_prism_face(local_boundary[i], order[0], order[1], gid, mesh, uv, points, weights);
797 normal_for_prism_face(local_boundary[i], normal);
798 break;
799 case BoundaryType::POLYGON:
800 quadrature_for_polygon_edge(local_boundary.element_id(), gid, order[0], mesh, uv, points, weights);
801 normal_for_polygon_edge(local_boundary.element_id(), gid, mesh, normal);
802 break;
803 case BoundaryType::INVALID:
804 assert(false);
805 break;
806 default:
807 assert(false);
808 }
809
810 normals.resize(points.rows(), normal.size());
811 for (int k = 0; k < normals.rows(); ++k)
812 normals.row(k) = normal;
813
814 return true;
815 }
816 } // namespace utils
817} // namespace polyfem
Quadrature quadrature
int edge_id
int x
Navigation::Index switch_edge(Navigation::Index idx) const override
Definition CMesh2D.hpp:87
int n_face_vertices(const int f_id) const override
number of vertices of a face
Definition CMesh2D.hpp:37
Navigation::Index get_index_from_face(int f, int lv=0) const override
Definition CMesh2D.hpp:83
virtual RowVectorNd point(const int global_index) const override
point coordinates
Definition CMesh2D.cpp:477
Boundary primitive IDs for a single element.
int element_id() const
Get the element's ID.
int size() const
Number of boundary primitives for the element.
int global_primitive_id(const int index) const
Get the i-th boundary primitive's global ID.
BoundaryType type() const
Get the type of boundary for the element.
Navigation::Index next_around_face(Navigation::Index idx) const
Definition Mesh2D.hpp:61
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:40
virtual double edge_length(const int gid) const
edge length
Definition Mesh.hpp:307
virtual double tri_area(const int gid) const
area of a tri face of a tet mesh
Definition Mesh.hpp:325
virtual double quad_area(const int gid) const
area of a quad face of an hex mesh
Definition Mesh.hpp:316
virtual bool is_volume() const =0
checks if mesh is volume
void get_quadrature(const int order, Quadrature &quad)
void get_quadrature(const int order, Quadrature &quad)
void get_quadrature(const int order, Quadrature &quad)
static void quadrature_for_polygon_edge(int face_id, int edge_id, int order, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
static Eigen::MatrixXd prism_local_node_coordinates_from_face(int lf)
static void quadrature_for_quad_face(int index, int order, const int gid, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
static bool sample_boundary(const mesh::LocalBoundary &local_boundary, const int n_samples, const mesh::Mesh &mesh, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples, Eigen::VectorXi &global_primitive_ids)
static void sample_parametric_prism_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void quadrature_for_quad_edge(int index, int order, const int gid, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
static void normal_for_quad_edge(int index, Eigen::MatrixXd &normal)
static void normal_for_tri_edge(int index, Eigen::MatrixXd &normal)
static void quadrature_for_tri_face(int index, int order, const int gid, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
static void normal_for_quad_face(int index, Eigen::MatrixXd &normal)
static void sample_parametric_tri_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static bool boundary_quadrature(const mesh::LocalBoundary &local_boundary, const QuadratureOrders &order, const mesh::Mesh &mesh, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::MatrixXd &normals, Eigen::VectorXd &weights, Eigen::VectorXi &global_primitive_ids)
static Eigen::MatrixXd hex_local_node_coordinates_from_face(int lf)
static void normal_for_prism_face(int index, Eigen::MatrixXd &normal)
static void normal_for_tri_face(int index, Eigen::MatrixXd &normal)
static void sample_parametric_quad_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void quadrature_for_prism_face(int index, int orderp, int orderq, const int gid, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
static void normal_for_polygon_edge(int face_id, int edge_id, const mesh::Mesh &mesh, Eigen::MatrixXd &normal)
static Eigen::MatrixXd tet_local_node_coordinates_from_face(int lf)
static void sample_parametric_quad_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void sample_polygon_edge(int face_id, int edge_id, int n_samples, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static void quadrature_for_tri_edge(int index, int order, const int gid, const mesh::Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::VectorXd &weights)
static void sample_parametric_tri_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
static Eigen::Matrix2d quad_local_node_coordinates_from_edge(int le)
static Eigen::Matrix2d tri_local_node_coordinates_from_edge(int le)
void q_nodes_2d(const int q, Eigen::MatrixXd &val)
void prism_nodes_3d(const int p, const int q, Eigen::MatrixXd &val)
void p_nodes_2d(const int p, Eigen::MatrixXd &val)
void p_nodes_3d(const int p, Eigen::MatrixXd &val)
void q_nodes_3d(const int q, Eigen::MatrixXd &val)
std::array< int, 2 > QuadratureOrders
Definition Types.hpp:19