19 #include "libmesh/libmesh_config.h"    21 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS    28 #include "libmesh/cell_inf_hex.h"    29 #include "libmesh/cell_inf_hex8.h"    30 #include "libmesh/face_quad4.h"    31 #include "libmesh/face_inf_quad4.h"    32 #include "libmesh/fe_type.h"    33 #include "libmesh/fe_interface.h"    34 #include "libmesh/inf_fe_map.h"    35 #include "libmesh/enum_elem_quality.h"    98   libmesh_assert_less (s, this->
n_sides());
   112   libmesh_assert_less (s, this->
n_sides());
   125                                      unsigned int side_node)
 const   127   libmesh_assert_less (side, this->
n_sides());
   136                                      unsigned int edge_node)
 const   138   libmesh_assert_less (edge, this->
n_edges());
   148   libmesh_assert_less (i, this->
n_sides());
   151   std::unique_ptr<Elem> face;
   169         face = std::make_unique<Quad4>();
   179         face = std::make_unique<InfQuad4>();
   184       libmesh_error_msg(
"Invalid side i = " << i);
   188   for (
auto n : face->node_index_range())
   197                        const unsigned int i)
   199   libmesh_assert_less (i, this->
n_sides());
   207         if (!side.get() || side->type() != 
QUAD4)
   221         if (!side.get() || side->type() != 
INFQUAD4)
   230       libmesh_error_msg(
"Invalid side i = " << i);
   236   for (
auto n : side->node_index_range())
   243                               const unsigned int s)
 const   246   libmesh_assert_less (s, this->
n_sides());
   248   return (s == 0 || c+1 == s || c == s%4);
   254                               const unsigned int s)
 const   256   libmesh_assert_less (e, this->
n_edges());
   257   libmesh_assert_less (s, this->
n_sides());
   266   libmesh_assert_less(e, this->
n_edges());
   280 std::vector<unsigned int>
   283   libmesh_assert_less(n, this->
n_nodes());
   292       auto trim = (n < 4) ? 0 : 2;
   325         const Real min = std::min(d02, d13);
   326         const Real max = std::max(d02, d13);
   328         libmesh_assert_not_equal_to (max, 0.0);
   353         std::vector<Real> edge_ratios(2);
   356         edge_ratios[0] = std::min(d01, d23) / std::max(d01, d23);
   357         edge_ratios[1] = std::min(d03, d12) / std::max(d03, d12);
   359         return *(std::min_element(edge_ratios.begin(), edge_ratios.end())) ;
   377         const Real sqrt3 = 1.73205080756888;
   384         const Real max_diag = std::max(d02, d13);
   386         libmesh_assert_not_equal_to ( max_diag, 0.0 );
   391         std::vector<Real> edges(4);
   392         edges[0]  = this->
length(0,1);
   393         edges[1]  = this->
length(1,2);
   394         edges[2]  = this->
length(2,3);
   395         edges[3]  = this->
length(0,3);
   397         const Real min_edge = *(std::min_element(edges.begin(), edges.end()));
   398         return sqrt3 * min_edge / max_diag ;
   415   libmesh_not_implemented();
   442     99,99,99,99,99,99,99,99, 
   452     99,99,99,99,99,99,99,99, 
   480   const Real tmp_min_distance_sq = std::min(pt0_o.
norm_sq(),
   489   const Real min_distance_sq = tmp_min_distance_sq
   495   const Real conservative_p_dist_sq = 1.01 * (
Point(p - my_origin).
norm_sq());
   499   if (conservative_p_dist_sq < min_distance_sq)
   507   Point p_o(p - my_origin);
   508   pt0_o /= pt0_o.
norm();
   509   pt1_o /= pt1_o.
norm();
   510   pt2_o /= pt2_o.
norm();
   511   pt3_o /= pt3_o.
norm();
   518                         (pt1_o - pt3_o).
norm_sq())*1.01;
   520   if ((p_o - pt0_o).
norm_sq() > max_h ||
   521       (p_o - pt1_o).
norm_sq() > max_h ||
   522       (p_o - pt2_o).
norm_sq() > max_h ||
   523       (p_o - pt3_o).
norm_sq() > max_h )
   542                                   const Real eps)
 const   544   const Real & xi = p(0);
   545   const Real & eta = p(1);
   546   const Real & zeta = p(2);
   549   return ((xi   >= -1.-eps) &&
   564 #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
static const int num_children
virtual std::vector< unsigned int > sides_on_edge(const unsigned int e) const override final
static const int num_sides
Geometric constants for all InfHexes. 
virtual dof_id_type low_order_key(const unsigned int s) const override
auto norm() const -> decltype(std::norm(T()))
virtual dof_id_type key() const
virtual bool is_face(const unsigned int i) const override final
We number faces last. 
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
static const int nodes_per_side
virtual bool is_edge_on_side(const unsigned int e, const unsigned int s) const override final
The libMesh namespace provides an interface to certain functionality in the library. 
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
static const unsigned short int _second_order_adjacent_vertices[8][2]
For higher-order elements, namely InfHex16 and InfHex18, the matrices for adjacent vertices of second...
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
static const int num_edges
static const Real _master_points[18][3]
Master element node locations. 
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
auto norm_sq() const -> decltype(std::norm(T()))
Real length(const unsigned int n1, const unsigned int n2) const
virtual unsigned int n_nodes() const =0
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
virtual bool is_edge(const unsigned int i) const override final
We number edges next. 
virtual unsigned int n_vertices() const override final
ElemQuality
Defines an enum for element quality metrics. 
virtual unsigned int n_children() const override final
static const unsigned int adjacent_edges_map[8][3]
This maps the  node to the 3 (or fewer) edge ids adjacent to the node. 
virtual bool is_vertex(const unsigned int i) const override final
We number vertices first. 
virtual bool is_flipped() const override final
virtual Point origin() const override
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the  node of the  side to element node numbers. 
static const unsigned int edge_sides_map[8][2]
This maps each edge to the sides that contain said edge. 
static const unsigned short int _second_order_vertex_child_index[18]
Vector that names the child vertex index for each second order node. 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
subdomain_id_type subdomain_id() const
virtual std::vector< unsigned int > edges_adjacent_to_node(const unsigned int n) const override
virtual bool contains_point(const Point &p, Real tol=TOLERANCE) const override
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override final
virtual bool on_reference_element(const Point &p, const Real eps=TOLERANCE) const override final
virtual Real quality(const ElemQuality q) const
virtual unsigned short dim() const override
virtual Real quality(const ElemQuality q) const override
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
This maps the  node of the  side to element node numbers. 
virtual unsigned int n_edges() const override final
static dof_id_type compute_key(dof_id_type n0)
virtual Order default_order() const =0
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
static const unsigned short int _second_order_vertex_child_number[18]
Vector that names a child sharing each second order node. 
static const int nodes_per_edge
virtual unsigned int n_sides() const override final
bool is_internal(const unsigned int i) const