23 Eigen::RowVector2d quad_local_node_coordinates(
int local_index)
25 assert(local_index >= 0 && local_index < 4);
28 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
31 Eigen::RowVector2d tri_local_node_coordinates(
int local_index)
35 return Eigen::RowVector2d(p(local_index, 0), p(local_index, 1));
38 Eigen::RowVector3d linear_tet_local_node_coordinates(
int local_index)
42 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
45 Eigen::RowVector3d linear_hex_local_node_coordinates(
int local_index)
49 return Eigen::RowVector3d(p(local_index, 0), p(local_index, 1), p(local_index, 2));
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);
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);
74 Eigen::Matrix<int, 4, 3> fv;
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));
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;
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));
106 auto endpoints = quad_local_node_coordinates_from_edge(index);
112 points.resize(quad.
points.rows(), endpoints.cols());
113 uv.resize(quad.
points.rows(), 2);
115 uv.col(0) = (1.0 - quad.
points.array());
116 uv.col(1) = quad.
points.array();
118 for (
int c = 0; c < 2; ++c)
120 points.col(c) = (1.0 - quad.
points.array()) * endpoints(0, c) + quad.
points.array() * endpoints(1, c);
129 auto endpoints = tri_local_node_coordinates_from_edge(index);
135 points.resize(quad.
points.rows(), endpoints.cols());
136 uv.resize(quad.
points.rows(), 2);
138 uv.col(0) = (1.0 - quad.
points.array());
139 uv.col(1) = quad.
points.array();
141 for (
int c = 0; c < 2; ++c)
143 points.col(c) = (1.0 - quad.
points.array()) * endpoints(0, c) + quad.
points.array() * endpoints(1, c);
152 auto endpoints = hex_local_node_coordinates_from_face(index);
158 const int n_pts = quad.
points.rows();
159 points.resize(n_pts, endpoints.cols());
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();
167 for (
int i = 0; i < n_pts; ++i)
169 const double b1 = quad.
points(i, 0);
170 const double b2 = 1 - b1;
172 const double b3 = quad.
points(i, 1);
173 const double b4 = 1 - b3;
175 for (
int c = 0; c < 3; ++c)
177 points(i, c) = b3 * (b1 * endpoints(0, c) + b2 * endpoints(1, c)) + b4 * (b1 * endpoints(3, c) + b2 * endpoints(2, c));
187 auto endpoints = tet_local_node_coordinates_from_face(index);
192 const int n_pts = quad.
points.rows();
193 points.resize(n_pts, endpoints.cols());
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();
200 for (
int i = 0; i < n_pts; ++i)
202 const double b1 = quad.
points(i, 0);
203 const double b3 = quad.
points(i, 1);
204 const double b2 = 1 - b1 - b3;
206 for (
int c = 0; c < 3; ++c)
208 points(i, c) = b1 * endpoints(0, c) + b2 * endpoints(1, c) + b3 * endpoints(2, c);
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());
223 uv.resize(n_samples, 2);
225 uv.col(0) = (1.0 - t.array());
226 uv.col(1) = t.array();
228 for (
int c = 0; c < 2; ++c)
230 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
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());
240 uv.resize(n_samples, 2);
242 uv.col(0) = (1.0 - t.array());
243 uv.col(1) = t.array();
245 for (
int c = 0; c < 2; ++c)
247 samples.col(c) = (1.0 - t.array()).matrix() * endpoints(0, c) + t * endpoints(1, c);
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());
259 uv.resize(n_samples * n_samples, endpoints.cols());
262 for (
int c = 0; c < 3; ++c)
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);
267 for (
int c = 0; c < 3; ++c)
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);
272 samples.col(c) = Eigen::Map<const Eigen::VectorXd>(
x.data(),
x.size());
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());
282 uv.resize(n_samples * n_samples, 3);
285 for (
int u = 0; u < n_samples; ++u)
287 for (
int v = 0; v < n_samples; ++v)
292 uv(counter, 0) = t(u);
293 uv(counter, 1) = t(v);
294 uv(counter, 2) = 1 - t(u) - t(v);
296 for (
int c = 0; c < 3; ++c)
298 samples(counter, c) = t(u) * endpoints(0, c) + t(v) * endpoints(1, c) + (1 - t(u) - t(v)) * endpoints(2, c);
303 samples.conservativeResize(counter, 3);
304 uv.conservativeResize(counter, 3);
329 auto p0 = mesh2d.
point(index.vertex);
331 const Eigen::VectorXd t = Eigen::VectorXd::LinSpaced(n_samples, 0, 1);
332 samples.resize(n_samples, p0.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()) * p0(c) + t.array() * p1(c);
367 auto p0 = mesh2d.
point(index.vertex);
373 points.resize(quad.
points.rows(), p0.cols());
374 uv.resize(quad.
points.rows(), 2);
376 uv.col(0) = (1.0 - quad.
points.array());
377 uv.col(1) = quad.
points.array();
379 for (
int c = 0; c < 2; ++c)
381 points.col(c) = (1.0 - quad.
points.array()) * p0(c) + quad.
points.array() * p1(c);
390 auto endpoints = quad_local_node_coordinates_from_edge(index);
391 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
400 auto endpoints = tri_local_node_coordinates_from_edge(index);
401 const Eigen::MatrixXd e = endpoints.row(0) - endpoints.row(1);
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);
416 if (index == 0 || index == 3 || index == 4)
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);
454 auto p0 = mesh2d.
point(index.vertex);
456 const Eigen::MatrixXd e = p0 - p1;
467 normals.resize(0, 0);
469 global_primitive_ids.resize(0);
471 for (
int i = 0; i < local_boundary.
size(); ++i)
474 Eigen::MatrixXd tmp_p, tmp_uv, tmp_n;
475 Eigen::VectorXd tmp_w;
476 switch (local_boundary.
type())
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);
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);
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);
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);
494 case BoundaryType::POLYGON:
498 case BoundaryType::INVALID:
505 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
506 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
508 points.conservativeResize(points.rows() + tmp_p.rows(), tmp_p.cols());
509 points.bottomRows(tmp_p.rows()) = tmp_p;
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;
515 weights.conservativeResize(weights.rows() + tmp_w.rows(), tmp_w.cols());
516 weights.bottomRows(tmp_w.rows()) = tmp_w;
518 global_primitive_ids.conservativeResize(global_primitive_ids.rows() + tmp_p.rows());
519 global_primitive_ids.bottomRows(tmp_p.rows()).setConstant(gid);
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());
533 samples.resize(0, 0);
534 global_primitive_ids.resize(0);
536 for (
int i = 0; i < local_boundary.
size(); ++i)
538 Eigen::MatrixXd tmp, tmp_uv;
539 switch (local_boundary.
type())
541 case BoundaryType::TRI_LINE:
542 sample_parametric_tri_edge(local_boundary[i], n_samples, tmp_uv, tmp);
544 case BoundaryType::QUAD_LINE:
545 sample_parametric_quad_edge(local_boundary[i], n_samples, tmp_uv, tmp);
547 case BoundaryType::QUAD:
548 sample_parametric_quad_face(local_boundary[i], n_samples, tmp_uv, tmp);
550 case BoundaryType::TRI:
551 sample_parametric_tri_face(local_boundary[i], n_samples, tmp_uv, tmp);
553 case BoundaryType::POLYGON:
556 case BoundaryType::INVALID:
563 uv.conservativeResize(uv.rows() + tmp_uv.rows(), tmp_uv.cols());
564 uv.bottomRows(tmp_uv.rows()) = tmp_uv;
566 samples.conservativeResize(samples.rows() + tmp.rows(), tmp.cols());
567 samples.bottomRows(tmp.rows()) = tmp;
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));
573 assert(uv.rows() == global_primitive_ids.size());
574 assert(samples.rows() == global_primitive_ids.size());
581 assert(local_boundary.
size() > i);
588 Eigen::MatrixXd normal;
589 switch (local_boundary.
type())
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);
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);
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);
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);
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);
611 case BoundaryType::INVALID:
618 normals.resize(points.rows(), normal.size());
619 for (
int k = 0; k < normals.rows(); ++k)
620 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 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)