Go to the documentation of this file.
   20 #ifndef LIBMESH_CELL_TET4_H 
   21 #define LIBMESH_CELL_TET4_H 
   24 #include "libmesh/cell_tet.h" 
   69   virtual ~Tet4() = 
default;
 
   84   virtual unsigned int n_sub_elem()
 const override { 
return 1; }
 
   89   virtual bool is_vertex(
const unsigned int i) 
const override;
 
   94   virtual bool is_edge(
const unsigned int i) 
const override;
 
   99   virtual bool is_face(
const unsigned int i) 
const override;
 
  106                                const unsigned int s) 
const override;
 
  108   virtual std::vector<unsigned int> 
nodes_on_side(
const unsigned int s) 
const override;
 
  115                                const unsigned int e) 
const override;
 
  122                                 const unsigned int s) 
const override;
 
  134   virtual bool is_linear ()
 const override { 
return true; }
 
  145   virtual std::unique_ptr<Elem> 
build_side_ptr (
const unsigned int i,
 
  146                                                 bool proxy=
true) 
override;
 
  152                                const unsigned int i) 
override;
 
  158   virtual std::unique_ptr<Elem> 
build_edge_ptr (
const unsigned int i) 
override;
 
  162                             std::vector<dof_id_type> & conn) 
const override;
 
  230 #ifdef LIBMESH_ENABLE_AMR 
  236                                   const unsigned int j,
 
  237                                   const unsigned int k) 
const override;
 
  254 #endif // LIBMESH_CELL_TET4_H 
  
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) 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.
 
virtual bool contains_point(const Point &p, Real tol) const override
Uses simple geometric tests to determine if the point p is inside the tetrahedron.
 
IOPackage
libMesh interfaces with several different software packages for the purposes of creating,...
 
The libMesh namespace provides an interface to certain functionality in the library.
 
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=true) override
Builds a TRI3 built coincident with face i.
 
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
 
virtual float embedding_matrix(const unsigned int i, const unsigned int j, const unsigned int k) const override
Matrix used to create the elements children.
 
static const int nodes_per_side
 
virtual unsigned int n_sub_elem() const override
 
virtual bool is_edge(const unsigned int i) const override
 
static const int nodes_per_edge
 
static const int num_children
 
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE2 built coincident with face i.
 
virtual dof_id_type key() const override
 
virtual bool is_linear() const override
 
Node * _nodelinks_data[num_nodes]
Data for links to nodes.
 
A Point defines a location in LIBMESH_DIM dimensional Real space.
 
virtual dof_id_type key() const
 
virtual unsigned int n_nodes() const override
 
A Node is like a Point, but with more information.
 
The Tet4 is an element in 3D composed of 4 nodes.
 
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
 
virtual Order default_order() const override
 
static const int num_edges
 
virtual Real volume() const override
An optimized method for computing the area of a 4-node tetrahedron.
 
static const float _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
 
static const int num_nodes
Geometric constants for Tet4.
 
virtual bool has_affine_map() const override
 
virtual bool is_face(const unsigned int i) const override
 
Tet4 & operator=(const Tet4 &)=delete
 
This is the base class from which all geometric element types are derived.
 
std::pair< Real, Real > min_and_max_angle() const
 
The Tet is an element in 3D composed of 4 sides.
 
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
 
virtual bool is_vertex(const unsigned int i) const override
 
LIBMESH_ENABLE_TOPOLOGY_CACHES
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
static const int num_sides
 
virtual ElemType type() const override
 
Tet4(Elem *p=nullptr)
Constructor.
 
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
This maps the  node of the  edge to element node numbers.
 
ElemType
Defines an enum for geometric element types.
 
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override