24 Eigen::RowVector2d quad_local_node_coordinates(
int local_index)
26 assert(local_index >= 0 && local_index < 4);
29 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
32 Eigen::RowVector2d tri_local_node_coordinates(
int local_index)
36 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
39 Eigen::RowVector3d linear_tet_local_node_coordinates(
int local_index)
43 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
46 Eigen::RowVector3d linear_hex_local_node_coordinates(
int local_index)
50 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
53 Eigen::RowVector3d linear_prism_local_node_coordinates(
int local_index)
57 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
60 Eigen::Matrix2d quad_local_node_coordinates_from_edge(
int le)
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);
69 Eigen::Matrix2d tri_local_node_coordinates_from_edge(
int le)
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);
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);
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);
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;
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));
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;
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));
131 Eigen::Matrix<int, 5, 4> fv;
132 fv.row(0) << 0, 1, 2, -1;
133 fv.row(1) << 3, 4, 5, -1;
135 fv.row(2) << 0, 1, 4, 3;
136 fv.row(3) << 1, 2, 5, 4;
137 fv.row(4) << 2, 0, 3, 5;
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));
149 auto endpoints = quad_local_node_coordinates_from_edge(index);
155 points.resize(quad.
points.rows(), endpoints.cols());
156 uv.resize(quad.
points.rows(), 2);
158 uv.col(0) = (1.0 - quad.
points.array());
159 uv.col(1) = quad.
points.array();
161 for (
int c = 0; c < 2; ++c)
163 points.col(c) = (1.0 - quad.
points.array()) * endpoints(0, c) + quad.
points.array() * endpoints(1, c);
172 auto endpoints = tri_local_node_coordinates_from_edge(index);
178 points.resize(quad.
points.rows(), endpoints.cols());
179 uv.resize(quad.
points.rows(), 2);
181 uv.col(0) = (1.0 - quad.
points.array());
182 uv.col(1) = quad.
points.array();
184 for (
int c = 0; c < 2; ++c)
186 points.col(c) = (1.0 - quad.
points.array()) * endpoints(0, c) + quad.
points.array() * endpoints(1, c);
195 auto endpoints = hex_local_node_coordinates_from_face(index);
201 const int n_pts = quad.
points.rows();
202 points.resize(n_pts, endpoints.cols());
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();
210 for (
int i = 0; i < n_pts; ++i)
212 const double b1 = quad.
points(i, 0);
213 const double b2 = 1 - b1;
215 const double b3 = quad.
points(i, 1);
216 const double b4 = 1 - b3;
218 for (
int c = 0; c < 3; ++c)
220 points(i, c) = b3 * (b1 * endpoints(0, c) + b2 * endpoints(1, c)) + b4 * (b1 * endpoints(3, c) + b2 * endpoints(2, c));
230 auto endpoints = tet_local_node_coordinates_from_face(index);
235 const int n_pts = quad.
points.rows();
236 points.resize(n_pts, endpoints.cols());
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();
243 for (
int i = 0; i < n_pts; ++i)
245 const double b1 = quad.
points(i, 0);
246 const double b3 = quad.
points(i, 1);
247 const double b2 = 1 - b1 - b3;
249 for (
int c = 0; c < 3; ++c)
251 points(i, c) = b1 * endpoints(0, c) + b2 * endpoints(1, c) + b3 * endpoints(2, c);
262 auto endpoints = prism_local_node_coordinates_from_face(index);
270 const int n_pts = quad.
points.rows();
271 points.resize(n_pts, endpoints.cols());
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();
278 for (
int i = 0; i < n_pts; ++i)
280 const double b1 = quad.
points(i, 0);
281 const double b3 = quad.
points(i, 1);
282 const double b2 = 1 - b1 - b3;
284 for (
int c = 0; c < 3; ++c)
286 points(i, c) = b1 * endpoints(0, c) + b2 * endpoints(1, c) + b3 * endpoints(2, c);
300 const int n_pts = quad.
points.rows();
301 points.resize(n_pts, endpoints.cols());
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();
309 for (
int i = 0; i < n_pts; ++i)
311 const double b1 = quad.
points(i, 0);
312 const double b2 = 1 - b1;
314 const double b3 = quad.
points(i, 1);
315 const double b4 = 1 - b3;
317 for (
int c = 0; c < 3; ++c)
319 points(i, c) = b3 * (b1 * endpoints(0, c) + b2 * endpoints(1, c)) + b4 * (b1 * endpoints(3, c) + b2 * endpoints(2, c));
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());
334 uv.resize(n_samples, 2);
336 uv.col(0) = (1.0 - t.array());
337 uv.col(1) = t.array();
339 for (
int c = 0; c < 2; ++c)
341 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
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());
351 uv.resize(n_samples, 2);
353 uv.col(0) = (1.0 - t.array());
354 uv.col(1) = t.array();
356 for (
int c = 0; c < 2; ++c)
358 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
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());
370 uv.resize(n_samples * n_samples, endpoints.cols());
373 for (
int c = 0; c < 3; ++c)
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);
378 for (
int c = 0; c < 3; ++c)
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);
383 samples.col(c) = Eigen::Map<const Eigen::VectorXd>(
x.data(),
x.size());
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());
393 uv.resize(n_samples * n_samples, 3);
396 for (
int u = 0; u < n_samples; ++u)
398 for (
int v = 0; v < n_samples; ++v)
403 uv(counter, 0) = t(u);
404 uv(counter, 1) = t(v);
405 uv(counter, 2) = 1 - t(u) - t(v);
407 for (
int c = 0; c < 3; ++c)
409 samples(counter, c) = t(u) * endpoints(0, c) + t(v) * endpoints(1, c) + (1 - t(u) - t(v)) * endpoints(2, c);
414 samples.conservativeResize(counter, 3);
415 uv.conservativeResize(counter, 3);
420 auto endpoints = prism_local_node_coordinates_from_face(index);
424 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
425 samples.resize(n_samples * n_samples, endpoints.cols());
427 uv.resize(n_samples * n_samples, 3);
430 for (
int u = 0; u < n_samples; ++u)
432 for (
int v = 0; v < n_samples; ++v)
437 uv(counter, 0) = t(u);
438 uv(counter, 1) = t(v);
439 uv(counter, 2) = 1 - t(u) - t(v);
441 for (
int c = 0; c < 3; ++c)
443 samples(counter, c) = t(u) * endpoints(0, c) + t(v) * endpoints(1, c) + (1 - t(u) - t(v)) * endpoints(2, c);
448 samples.conservativeResize(counter, 3);
449 uv.conservativeResize(counter, 3);
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());
458 uv.resize(n_samples * n_samples, endpoints.cols());
461 for (
int c = 0; c < 3; ++c)
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);
466 for (
int c = 0; c < 3; ++c)
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);
471 samples.col(c) = Eigen::Map<const Eigen::VectorXd>(
x.data(),
x.size());
498 auto p0 = mesh2d.
point(index.vertex);
500 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
501 samples.resize(n_samples, p0.cols());
503 uv.resize(n_samples, 2);
505 uv.col(0) = (1.0 - t.array());
506 uv.col(1) = t.array();
508 for (
int c = 0; c < 2; ++c)
510 samples.col(c) = (1.0 - t.array()) * p0(c) + t.array() * p1(c);
536 auto p0 = mesh2d.
point(index.vertex);
542 points.resize(quad.
points.rows(), p0.cols());
543 uv.resize(quad.
points.rows(), 2);
545 uv.col(0) = (1.0 - quad.
points.array());
546 uv.col(1) = quad.
points.array();
548 for (
int c = 0; c < 2; ++c)
550 points.col(c) = (1.0 - quad.
points.array()) * p0(c) + quad.
points.array() * p1(c);
559 auto endpoints = quad_local_node_coordinates_from_edge(index);
560 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
569 auto endpoints = tri_local_node_coordinates_from_edge(index);
570 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
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);
585 if (index == 0 || index == 3 || index == 4)
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);
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);
635 auto p0 = mesh2d.
point(index.vertex);
637 const Eigen::MatrixXd e = p0 - p1;
648 normals.resize(0, 0);
650 global_primitive_ids.resize(0);
652 for (
int i = 0; i < local_boundary.
size(); ++i)
655 Eigen::MatrixXd tmp_p, tmp_uv, tmp_n;
656 Eigen::VectorXd tmp_w;
657 switch (local_boundary.
type())
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);
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);
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);
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);
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);
679 case BoundaryType::POLYGON:
683 case BoundaryType::INVALID:
690 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
691 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
693 points.conservativeResize(points.rows() + tmp_p.rows(), tmp_p.cols());
694 points.bottomRows(tmp_p.rows()) = tmp_p;
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;
700 weights.conservativeResize(weights.rows() + tmp_w.rows(), tmp_w.cols());
701 weights.bottomRows(tmp_w.rows()) = tmp_w;
703 global_primitive_ids.conservativeResize(global_primitive_ids.rows() + tmp_p.rows());
704 global_primitive_ids.bottomRows(tmp_p.rows()).setConstant(gid);
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());
718 samples.resize(0, 0);
719 global_primitive_ids.resize(0);
721 for (
int i = 0; i < local_boundary.
size(); ++i)
723 Eigen::MatrixXd tmp, tmp_uv;
724 switch (local_boundary.
type())
726 case BoundaryType::TRI_LINE:
727 sample_parametric_tri_edge(local_boundary[i], n_samples, tmp_uv, tmp);
729 case BoundaryType::QUAD_LINE:
730 sample_parametric_quad_edge(local_boundary[i], n_samples, tmp_uv, tmp);
732 case BoundaryType::QUAD:
733 sample_parametric_quad_face(local_boundary[i], n_samples, tmp_uv, tmp);
735 case BoundaryType::TRI:
736 sample_parametric_tri_face(local_boundary[i], n_samples, tmp_uv, tmp);
738 case BoundaryType::PRISM:
739 sample_parametric_prism_face(local_boundary[i], n_samples, tmp_uv, tmp);
741 case BoundaryType::POLYGON:
744 case BoundaryType::INVALID:
751 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
752 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
754 samples.conservativeResize(samples.rows() + tmp.rows(), tmp.cols());
755 samples.bottomRows(tmp.rows()) = tmp;
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));
761 assert(uv.rows() == global_primitive_ids.size());
762 assert(samples.rows() == global_primitive_ids.size());
769 assert(local_boundary.
size() > i);
776 Eigen::MatrixXd normal;
777 switch (local_boundary.
type())
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);
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);
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);
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);
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);
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);
803 case BoundaryType::INVALID:
810 normals.resize(points.rows(), normal.size());
811 for (
int k = 0; k < normals.rows(); ++k)
812 normals.row(k) = normal;
Navigation::Index switch_edge(Navigation::Index idx) const override
int n_face_vertices(const int f_id) const override
number of vertices of a face
Navigation::Index get_index_from_face(int f, int lv=0) const override
virtual RowVectorNd point(const int global_index) const override
point coordinates
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
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
virtual double edge_length(const int gid) const
edge length
virtual double tri_area(const int gid) const
area of a tri face of a tet mesh
virtual double quad_area(const int gid) const
area of a quad face of an hex mesh
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