PolyFEM
Loading...
Searching...
No Matches
BoundarySampler.cpp
Go to the documentation of this file.
1#include "BoundarySampler.hpp"
2
6
9
12
13#include <cassert>
14
15namespace polyfem
16{
17 using namespace mesh;
18 using namespace quadrature;
19 namespace utils
20 {
21 namespace
22 {
23 Eigen::RowVector2d quad_local_node_coordinates(int local_index)
24 {
25 assert(local_index >= 0 && local_index < 4);
26 Eigen::MatrixXd p;
28 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
29 }
30
31 Eigen::RowVector2d tri_local_node_coordinates(int local_index)
32 {
33 Eigen::MatrixXd p;
35 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
36 }
37
38 Eigen::RowVector3d linear_tet_local_node_coordinates(int local_index)
39 {
40 Eigen::MatrixXd p;
42 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
43 }
44
45 Eigen::RowVector3d linear_hex_local_node_coordinates(int local_index)
46 {
47 Eigen::MatrixXd p;
49 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
50 }
51
52 } // namespace
53
55 {
56 Eigen::Matrix2d res(2, 2);
57 res.row(0) = quad_local_node_coordinates(le);
58 res.row(1) = quad_local_node_coordinates((le + 1) % 4);
59
60 return res;
61 }
62
64 {
65 Eigen::Matrix2d res(2, 2);
66 res.row(0) = tri_local_node_coordinates(le);
67 res.row(1) = tri_local_node_coordinates((le + 1) % 3);
68
69 return res;
70 }
71
73 {
74 Eigen::Matrix<int, 4, 3> fv;
75 fv.row(0) << 0, 1, 2;
76 fv.row(1) << 0, 1, 3;
77 fv.row(2) << 1, 2, 3;
78 fv.row(3) << 2, 0, 3;
79
80 Eigen::MatrixXd res(3, 3);
81 for (int i = 0; i < 3; ++i)
82 res.row(i) = linear_tet_local_node_coordinates(fv(lf, i));
83
84 return res;
85 }
86
88 {
89 Eigen::Matrix<int, 6, 4> fv;
90 fv.row(0) << 0, 3, 7, 4;
91 fv.row(1) << 1, 2, 6, 5;
92 fv.row(2) << 0, 1, 5, 4;
93 fv.row(3) << 3, 2, 6, 7;
94 fv.row(4) << 0, 1, 2, 3;
95 fv.row(5) << 4, 5, 6, 7;
96
97 Eigen::MatrixXd res(4, 3);
98 for (int i = 0; i < 4; ++i)
99 res.row(i) = linear_hex_local_node_coordinates(fv(lf, i));
100
101 return res;
102 }
103
104 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)
105 {
106 auto endpoints = quad_local_node_coordinates_from_edge(index);
107
108 Quadrature quad;
109 LineQuadrature quad_rule;
110 quad_rule.get_quadrature(order, quad);
111
112 points.resize(quad.points.rows(), endpoints.cols());
113 uv.resize(quad.points.rows(), 2);
114
115 uv.col(0) = (1.0 - quad.points.array());
116 uv.col(1) = quad.points.array();
117
118 for (int c = 0; c < 2; ++c)
119 {
120 points.col(c) = (1.0 - quad.points.array()) * endpoints(0, c) + quad.points.array() * endpoints(1, c);
121 }
122
123 weights = quad.weights;
124 weights *= mesh.edge_length(gid);
125 }
126
127 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)
128 {
129 auto endpoints = tri_local_node_coordinates_from_edge(index);
130
131 Quadrature quad;
132 LineQuadrature quad_rule;
133 quad_rule.get_quadrature(order, quad);
134
135 points.resize(quad.points.rows(), endpoints.cols());
136 uv.resize(quad.points.rows(), 2);
137
138 uv.col(0) = (1.0 - quad.points.array());
139 uv.col(1) = quad.points.array();
140
141 for (int c = 0; c < 2; ++c)
142 {
143 points.col(c) = (1.0 - quad.points.array()) * endpoints(0, c) + quad.points.array() * endpoints(1, c);
144 }
145
146 weights = quad.weights;
147 weights *= mesh.edge_length(gid);
148 }
149
150 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)
151 {
152 auto endpoints = hex_local_node_coordinates_from_face(index);
153
154 Quadrature quad;
155 QuadQuadrature quad_rule;
156 quad_rule.get_quadrature(order, quad);
157
158 const int n_pts = quad.points.rows();
159 points.resize(n_pts, endpoints.cols());
160
161 uv.resize(quad.points.rows(), 4);
162 uv.col(0) = quad.points.col(0);
163 uv.col(1) = 1 - uv.col(0).array();
164 uv.col(2) = quad.points.col(1);
165 uv.col(3) = 1 - uv.col(2).array();
166
167 for (int i = 0; i < n_pts; ++i)
168 {
169 const double b1 = quad.points(i, 0);
170 const double b2 = 1 - b1;
171
172 const double b3 = quad.points(i, 1);
173 const double b4 = 1 - b3;
174
175 for (int c = 0; c < 3; ++c)
176 {
177 points(i, c) = b3 * (b1 * endpoints(0, c) + b2 * endpoints(1, c)) + b4 * (b1 * endpoints(3, c) + b2 * endpoints(2, c));
178 }
179 }
180
181 weights = quad.weights;
182 weights *= mesh.quad_area(gid);
183 }
184
185 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)
186 {
187 auto endpoints = tet_local_node_coordinates_from_face(index);
188 Quadrature quad;
189 TriQuadrature quad_rule;
190 quad_rule.get_quadrature(order, quad);
191
192 const int n_pts = quad.points.rows();
193 points.resize(n_pts, endpoints.cols());
194
195 uv.resize(quad.points.rows(), 3);
196 uv.col(0) = quad.points.col(0);
197 uv.col(1) = quad.points.col(1);
198 uv.col(2) = 1 - uv.col(0).array() - uv.col(1).array();
199
200 for (int i = 0; i < n_pts; ++i)
201 {
202 const double b1 = quad.points(i, 0);
203 const double b3 = quad.points(i, 1);
204 const double b2 = 1 - b1 - b3;
205
206 for (int c = 0; c < 3; ++c)
207 {
208 points(i, c) = b1 * endpoints(0, c) + b2 * endpoints(1, c) + b3 * endpoints(2, c);
209 }
210 }
211
212 weights = quad.weights;
213 // 2 * because weights sum to 1/2 already
214 weights *= 2 * mesh.tri_area(gid);
215 }
216
217 void utils::BoundarySampler::sample_parametric_quad_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
218 {
219 auto endpoints = quad_local_node_coordinates_from_edge(index);
220 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
221 samples.resize(n_samples, endpoints.cols());
222
223 uv.resize(n_samples, 2);
224
225 uv.col(0) = (1.0 - t.array());
226 uv.col(1) = t.array();
227
228 for (int c = 0; c < 2; ++c)
229 {
230 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
231 }
232 }
233
234 void utils::BoundarySampler::sample_parametric_tri_edge(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
235 {
236 auto endpoints = tri_local_node_coordinates_from_edge(index);
237 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
238 samples.resize(n_samples, endpoints.cols());
239
240 uv.resize(n_samples, 2);
241
242 uv.col(0) = (1.0 - t.array());
243 uv.col(1) = t.array();
244
245 for (int c = 0; c < 2; ++c)
246 {
247 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
248 }
249 }
250
251 void utils::BoundarySampler::sample_parametric_quad_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
252 {
253 auto endpoints = hex_local_node_coordinates_from_face(index);
254 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
255 samples.resize(n_samples * n_samples, endpoints.cols());
256 Eigen::MatrixXd left(n_samples, endpoints.cols());
257 Eigen::MatrixXd right(n_samples, endpoints.cols());
258
259 uv.resize(n_samples * n_samples, endpoints.cols());
260 uv.setZero();
261
262 for (int c = 0; c < 3; ++c)
263 {
264 left.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
265 right.col(c) = (1.0 - t.array()).matrix() * endpoints(3, c) + t * endpoints(2, c);
266 }
267 for (int c = 0; c < 3; ++c)
268 {
269 Eigen::MatrixXd x = (1.0 - t.array()).matrix() * left.col(c).transpose() + t * right.col(c).transpose();
270 assert(x.size() == n_samples * n_samples);
271
272 samples.col(c) = Eigen::Map<const Eigen::VectorXd>(x.data(), x.size());
273 }
274 }
275
276 void utils::BoundarySampler::sample_parametric_tri_face(int index, int n_samples, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
277 {
278 auto endpoints = tet_local_node_coordinates_from_face(index);
279 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
280 samples.resize(n_samples * n_samples, endpoints.cols());
281
282 uv.resize(n_samples * n_samples, 3);
283
284 int counter = 0;
285 for (int u = 0; u < n_samples; ++u)
286 {
287 for (int v = 0; v < n_samples; ++v)
288 {
289 if (t(u) + t(v) > 1)
290 continue;
291
292 uv(counter, 0) = t(u);
293 uv(counter, 1) = t(v);
294 uv(counter, 2) = 1 - t(u) - t(v);
295
296 for (int c = 0; c < 3; ++c)
297 {
298 samples(counter, c) = t(u) * endpoints(0, c) + t(v) * endpoints(1, c) + (1 - t(u) - t(v)) * endpoints(2, c);
299 }
300 ++counter;
301 }
302 }
303 samples.conservativeResize(counter, 3);
304 uv.conservativeResize(counter, 3);
305 }
306
307 void utils::BoundarySampler::sample_polygon_edge(int face_id, int edge_id, int n_samples, const Mesh &mesh, Eigen::MatrixXd &uv, Eigen::MatrixXd &samples)
308 {
309 assert(!mesh.is_volume());
310
311 const CMesh2D &mesh2d = dynamic_cast<const CMesh2D &>(mesh);
312
313 auto index = mesh2d.get_index_from_face(face_id);
314
315 bool found = false;
316 for (int i = 0; i < mesh2d.n_face_vertices(face_id); ++i)
317 {
318 if (index.edge == edge_id)
319 {
320 found = true;
321 break;
322 }
323
324 index = mesh2d.next_around_face(index);
325 }
326
327 assert(found);
328
329 auto p0 = mesh2d.point(index.vertex);
330 auto p1 = mesh2d.point(mesh2d.switch_edge(index).vertex);
331 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
332 samples.resize(n_samples, p0.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()) * p0(c) + t.array() * p1(c);
342 }
343 }
344
345 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)
346 {
347 assert(!mesh.is_volume());
348
349 const CMesh2D &mesh2d = dynamic_cast<const CMesh2D &>(mesh);
350
351 auto index = mesh2d.get_index_from_face(face_id);
352
353 bool found = false;
354 for (int i = 0; i < mesh2d.n_face_vertices(face_id); ++i)
355 {
356 if (index.edge == edge_id)
357 {
358 found = true;
359 break;
360 }
361
362 index = mesh2d.next_around_face(index);
363 }
364
365 assert(found);
366
367 auto p0 = mesh2d.point(index.vertex);
368 auto p1 = mesh2d.point(mesh2d.switch_edge(index).vertex);
369 Quadrature quad;
370 LineQuadrature quad_rule;
371 quad_rule.get_quadrature(order, quad);
372
373 points.resize(quad.points.rows(), p0.cols());
374 uv.resize(quad.points.rows(), 2);
375
376 uv.col(0) = (1.0 - quad.points.array());
377 uv.col(1) = quad.points.array();
378
379 for (int c = 0; c < 2; ++c)
380 {
381 points.col(c) = (1.0 - quad.points.array()) * p0(c) + quad.points.array() * p1(c);
382 }
383
384 weights = quad.weights;
385 weights *= mesh.edge_length(edge_id);
386 }
387
388 void utils::BoundarySampler::normal_for_quad_edge(int index, Eigen::MatrixXd &normal)
389 {
390 auto endpoints = quad_local_node_coordinates_from_edge(index);
391 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
392 normal.resize(1, 2);
393 normal(0) = -e(1);
394 normal(1) = e(0);
395 normal.normalize();
396 }
397
398 void utils::BoundarySampler::normal_for_tri_edge(int index, Eigen::MatrixXd &normal)
399 {
400 auto endpoints = tri_local_node_coordinates_from_edge(index);
401 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
402 normal.resize(1, 2);
403 normal(0) = -e(1);
404 normal(1) = e(0);
405 normal.normalize();
406 }
407
408 void utils::BoundarySampler::normal_for_quad_face(int index, Eigen::MatrixXd &normal)
409 {
410 auto endpoints = hex_local_node_coordinates_from_face(index);
411 const Eigen::RowVector3d e1 = endpoints.row(0) - endpoints.row(1);
412 const Eigen::RowVector3d e2 = endpoints.row(0) - endpoints.row(2);
413 normal = e1.cross(e2);
414 normal.normalize();
415
416 if (index == 0 || index == 3 || index == 4)
417 normal *= -1;
418 }
419
420 void utils::BoundarySampler::normal_for_tri_face(int index, Eigen::MatrixXd &normal)
421 {
422 auto endpoints = tet_local_node_coordinates_from_face(index);
423 const Eigen::RowVector3d e1 = endpoints.row(0) - endpoints.row(1);
424 const Eigen::RowVector3d e2 = endpoints.row(0) - endpoints.row(2);
425 normal = e1.cross(e2);
426 normal.normalize();
427
428 if (index == 0)
429 normal *= -1;
430 }
431
432 void utils::BoundarySampler::normal_for_polygon_edge(int face_id, int edge_id, const Mesh &mesh, Eigen::MatrixXd &normal)
433 {
434 assert(!mesh.is_volume());
435
436 const CMesh2D &mesh2d = dynamic_cast<const CMesh2D &>(mesh);
437
438 auto index = mesh2d.get_index_from_face(face_id);
439
440 bool found = false;
441 for (int i = 0; i < mesh2d.n_face_vertices(face_id); ++i)
442 {
443 if (index.edge == edge_id)
444 {
445 found = true;
446 break;
447 }
448
449 index = mesh2d.next_around_face(index);
450 }
451
452 assert(found);
453
454 auto p0 = mesh2d.point(index.vertex);
455 auto p1 = mesh2d.point(mesh2d.switch_edge(index).vertex);
456 const Eigen::MatrixXd e = p0 - p1;
457 normal.resize(1, 2);
458 normal(0) = -e(1);
459 normal(1) = e(0);
460 normal.normalize();
461 }
462
463 bool utils::BoundarySampler::boundary_quadrature(const LocalBoundary &local_boundary, const int 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)
464 {
465 uv.resize(0, 0);
466 points.resize(0, 0);
467 normals.resize(0, 0);
468 weights.resize(0);
469 global_primitive_ids.resize(0);
470
471 for (int i = 0; i < local_boundary.size(); ++i)
472 {
473 const int gid = local_boundary.global_primitive_id(i);
474 Eigen::MatrixXd tmp_p, tmp_uv, tmp_n;
475 Eigen::VectorXd tmp_w;
476 switch (local_boundary.type())
477 {
478 case BoundaryType::TRI_LINE:
479 quadrature_for_tri_edge(local_boundary[i], order, gid, mesh, tmp_uv, tmp_p, tmp_w);
480 normal_for_tri_edge(local_boundary[i], tmp_n);
481 break;
482 case BoundaryType::QUAD_LINE:
483 quadrature_for_quad_edge(local_boundary[i], order, gid, mesh, tmp_uv, tmp_p, tmp_w);
484 normal_for_quad_edge(local_boundary[i], tmp_n);
485 break;
486 case BoundaryType::QUAD:
487 quadrature_for_quad_face(local_boundary[i], order, gid, mesh, tmp_uv, tmp_p, tmp_w);
488 normal_for_quad_face(local_boundary[i], tmp_n);
489 break;
490 case BoundaryType::TRI:
491 quadrature_for_tri_face(local_boundary[i], order, gid, mesh, tmp_uv, tmp_p, tmp_w);
492 normal_for_tri_face(local_boundary[i], tmp_n);
493 break;
494 case BoundaryType::POLYGON:
495 quadrature_for_polygon_edge(local_boundary.element_id(), local_boundary.global_primitive_id(i), order, mesh, tmp_uv, tmp_p, tmp_w);
496 normal_for_polygon_edge(local_boundary.element_id(), local_boundary.global_primitive_id(i), mesh, tmp_n);
497 break;
498 case BoundaryType::INVALID:
499 assert(false);
500 break;
501 default:
502 assert(false);
503 }
504
505 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
506 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
507
508 points.conservativeResize(points.rows() + tmp_p.rows(), tmp_p.cols());
509 points.bottomRows(tmp_p.rows()) = tmp_p;
510
511 normals.conservativeResize(normals.rows() + tmp_p.rows(), tmp_p.cols());
512 for (int k = normals.rows() - tmp_p.rows(); k < normals.rows(); ++k)
513 normals.row(k) = tmp_n;
514
515 weights.conservativeResize(weights.rows() + tmp_w.rows(), tmp_w.cols());
516 weights.bottomRows(tmp_w.rows()) = tmp_w;
517
518 global_primitive_ids.conservativeResize(global_primitive_ids.rows() + tmp_p.rows());
519 global_primitive_ids.bottomRows(tmp_p.rows()).setConstant(gid);
520 }
521
522 assert(uv.rows() == global_primitive_ids.size());
523 assert(points.rows() == global_primitive_ids.size());
524 assert(normals.rows() == global_primitive_ids.size());
525 assert(weights.size() == global_primitive_ids.size());
526
527 return true;
528 }
529
530 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)
531 {
532 uv.resize(0, 0);
533 samples.resize(0, 0);
534 global_primitive_ids.resize(0);
535
536 for (int i = 0; i < local_boundary.size(); ++i)
537 {
538 Eigen::MatrixXd tmp, tmp_uv;
539 switch (local_boundary.type())
540 {
541 case BoundaryType::TRI_LINE:
542 sample_parametric_tri_edge(local_boundary[i], n_samples, tmp_uv, tmp);
543 break;
544 case BoundaryType::QUAD_LINE:
545 sample_parametric_quad_edge(local_boundary[i], n_samples, tmp_uv, tmp);
546 break;
547 case BoundaryType::QUAD:
548 sample_parametric_quad_face(local_boundary[i], n_samples, tmp_uv, tmp);
549 break;
550 case BoundaryType::TRI:
551 sample_parametric_tri_face(local_boundary[i], n_samples, tmp_uv, tmp);
552 break;
553 case BoundaryType::POLYGON:
554 sample_polygon_edge(local_boundary.element_id(), local_boundary.global_primitive_id(i), n_samples, mesh, tmp_uv, tmp);
555 break;
556 case BoundaryType::INVALID:
557 assert(false);
558 break;
559 default:
560 assert(false);
561 }
562
563 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
564 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
565
566 samples.conservativeResize(samples.rows() + tmp.rows(), tmp.cols());
567 samples.bottomRows(tmp.rows()) = tmp;
568
569 global_primitive_ids.conservativeResize(global_primitive_ids.rows() + tmp.rows());
570 global_primitive_ids.bottomRows(tmp.rows()).setConstant(local_boundary.global_primitive_id(i));
571 }
572
573 assert(uv.rows() == global_primitive_ids.size());
574 assert(samples.rows() == global_primitive_ids.size());
575
576 return true;
577 }
578
579 bool utils::BoundarySampler::boundary_quadrature(const mesh::LocalBoundary &local_boundary, const int order, const mesh::Mesh &mesh, const int i, const bool skip_computation, Eigen::MatrixXd &uv, Eigen::MatrixXd &points, Eigen::MatrixXd &normals, Eigen::VectorXd &weights)
580 {
581 assert(local_boundary.size() > i);
582
583 uv.resize(0, 0);
584 points.resize(0, 0);
585 weights.resize(0);
586 const int gid = local_boundary.global_primitive_id(i);
587
588 Eigen::MatrixXd normal;
589 switch (local_boundary.type())
590 {
591 case BoundaryType::TRI_LINE:
592 quadrature_for_tri_edge(local_boundary[i], order, gid, mesh, uv, points, weights);
593 normal_for_tri_edge(local_boundary[i], normal);
594 break;
595 case BoundaryType::QUAD_LINE:
596 quadrature_for_quad_edge(local_boundary[i], order, gid, mesh, uv, points, weights);
597 normal_for_quad_edge(local_boundary[i], normal);
598 break;
599 case BoundaryType::QUAD:
600 quadrature_for_quad_face(local_boundary[i], order, gid, mesh, uv, points, weights);
601 normal_for_quad_face(local_boundary[i], normal);
602 break;
603 case BoundaryType::TRI:
604 quadrature_for_tri_face(local_boundary[i], order, gid, mesh, uv, points, weights);
605 normal_for_tri_face(local_boundary[i], normal);
606 break;
607 case BoundaryType::POLYGON:
608 quadrature_for_polygon_edge(local_boundary.element_id(), gid, order, mesh, uv, points, weights);
609 normal_for_polygon_edge(local_boundary.element_id(), gid, mesh, normal);
610 break;
611 case BoundaryType::INVALID:
612 assert(false);
613 break;
614 default:
615 assert(false);
616 }
617
618 normals.resize(points.rows(), normal.size());
619 for (int k = 0; k < normals.rows(); ++k)
620 normals.row(k) = normal;
621
622 return true;
623 }
624 } // namespace utils
625} // 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:39
virtual double edge_length(const int gid) const
edge length
Definition Mesh.hpp:306
virtual double tri_area(const int gid) const
area of a tri face of a tet mesh
Definition Mesh.hpp:324
virtual double quad_area(const int gid) const
area of a quad face of an hex mesh
Definition Mesh.hpp:315
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 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 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 Eigen::MatrixXd hex_local_node_coordinates_from_face(int lf)
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 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 bool boundary_quadrature(const mesh::LocalBoundary &local_boundary, const int 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 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 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)