19 #include "libmesh/face_tri.h" 20 #include "libmesh/edge_edge2.h" 21 #include "libmesh/face_tri3.h" 22 #include "libmesh/enum_elem_quality.h" 61 libmesh_assert_less (s, this->
n_sides());
71 libmesh_assert_less (s, this->
n_sides());
80 unsigned int side_node)
const 82 libmesh_assert_less (side, this->
n_sides());
91 unsigned int edge_node)
const 109 libmesh_assert_less (i, this->
n_sides());
111 std::unique_ptr<Elem> edge = std::make_unique<Edge2>();
113 for (
auto n : edge->node_index_range())
122 const unsigned int i)
124 this->simple_side_ptr<Tri,Tri3>(side, i,
EDGE2);
130 const unsigned int s)
const 133 libmesh_assert_less (s, this->
n_sides());
135 return (c == s || c == (s+1)%3);
145 !this->
point(2)(2) &&
154 std::vector<unsigned int>
157 libmesh_assert_less(n, this->
n_nodes());
188 static const unsigned int opposite[3] = {1, 2, 0};
191 static const unsigned int other[3] = {2, 0, 1};
194 auto vertex_aspect_ratio = [&](
unsigned int i) ->
Real 200 Point v1 = m[other[i]] - m[i];
210 if (v0_norm == 0. || v1_norm == 0.)
215 Real v0s = v0_norm*sin_theta;
216 Real v1s = v1_norm*sin_theta;
220 auto [min0, max0] = std::minmax(v0_norm, v1s);
221 auto [min1, max1] = std::minmax(v0s, v1_norm);
224 return std::max(max0/min0, max1/min1);
227 return std::max(std::max(vertex_aspect_ratio(0), vertex_aspect_ratio(1)), vertex_aspect_ratio(2)) / std::sqrt(3);
248 if ((l1 <=0.) || (l2 <= 0.) || (l3 <= 0.))
251 const Real s1 = std::sin(std::acos(v1*v2/l1/l2)/2.);
253 const Real s2 = std::sin(std::acos(v1*v3/l1/l3)/2.);
254 const Real s3 = std::sin(std::acos(v2*v3/l2/l3)/2.);
256 return 8. * s1 * s2 * s3;
275 std::array<Real, 6> A =
285 lambda11 = A[0]*A[0] + A[2]*A[2] + A[4]*A[4],
286 lambda12 = A[0]*A[1] + A[2]*A[3] + A[4]*A[5],
287 lambda22 = A[1]*A[1] + A[3]*A[3] + A[5]*A[5];
291 Real den = lambda11 + lambda22 - lambda12;
296 Real alpha = std::sqrt(lambda11 * lambda22 - lambda12 * lambda12);
299 return std::sqrt(3) * alpha / den;
317 std::pair<Real, Real> bounds;
351 bounds.second = 1.155;
366 libMesh::out <<
"Warning: Invalid quality measure chosen." << std::endl;
376 const Real eps)
const 378 const Real & xi = p(0);
379 const Real & eta = p(1);
380 return ((xi >= 0.-eps) &&
382 ((xi + eta) <= 1.+eps));
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
static const Real _master_points[6][3]
Master element node locations.
auto norm() const -> decltype(std::norm(T()))
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
Calls local_side_node(edge, edge_node).
virtual bool is_face(const unsigned int i) const =0
static const int num_sides
Geometric constants for all Tris.
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Calls cross_norm_sq() and takes the square root of the result.
virtual bool is_flipped() const override final
virtual unsigned int n_vertices() const override final
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
virtual unsigned int n_sides() const override final
virtual std::vector< unsigned int > edges_adjacent_to_node(const unsigned int n) const override
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override final
virtual dof_id_type key() const override final
virtual dof_id_type low_order_key(const unsigned int s) const override
ElemQuality
Defines an enum for element quality metrics.
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
static const int nodes_per_side
virtual unsigned int n_nodes() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool on_reference_element(const Point &p, const Real eps=TOLERANCE) const override final
virtual bool is_vertex(const unsigned int i) const =0
virtual Real quality(const ElemQuality q) const
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
virtual unsigned int n_children() const override final
static dof_id_type compute_key(dof_id_type n0)
virtual Real quality(const ElemQuality q) const override
A Point defines a location in LIBMESH_DIM dimensional Real space.
dof_id_type node_id(const unsigned int i) const
const Point & point(const unsigned int i) const
virtual bool is_edge(const unsigned int i) const =0
static const unsigned int adjacent_sides_map[3][2]
This maps the node to the (in this case) 2 side ids adjacent to the node.