20 #ifndef LIBMESH_CELL_TET4_H 21 #define LIBMESH_CELL_TET4_H 24 #include "libmesh/cell_tet.h" 75 virtual ~Tet4() =
default;
90 virtual unsigned int n_sub_elem()
const override {
return 1; }
95 virtual bool is_vertex(
const unsigned int i)
const override;
100 virtual bool is_edge(
const unsigned int i)
const override;
105 virtual bool is_face(
const unsigned int i)
const override;
112 const unsigned int s)
const override;
114 virtual std::vector<unsigned int>
nodes_on_side(
const unsigned int s)
const override;
116 virtual std::vector<unsigned int>
nodes_on_edge(
const unsigned int e)
const override;
123 const unsigned int e)
const override;
130 const unsigned int s)
const override;
147 virtual bool is_linear ()
const override {
return true; }
158 virtual std::unique_ptr<Elem>
build_side_ptr (
const unsigned int i)
override;
164 const unsigned int i)
override;
173 virtual std::unique_ptr<Elem>
build_edge_ptr (
const unsigned int i)
override;
178 virtual void build_edge_ptr (std::unique_ptr<Elem> & edge,
const unsigned int i)
override;
182 std::vector<dof_id_type> & conn)
const override;
244 virtual void permute(
unsigned int perm_num)
override final;
259 #ifdef LIBMESH_ENABLE_AMR 265 const unsigned int j,
266 const unsigned int k)
const override;
283 #endif // LIBMESH_CELL_TET4_H virtual void permute(unsigned int perm_num) override final
Permutes the element (by swapping node and neighbor pointers) according to the specified index...
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
A Node is like a Point, but with more information.
virtual Point true_centroid() const override
The centroid of a 4-node tetrahedron is simply given by the average of its vertex positions...
std::pair< Real, Real > min_and_max_angle() const
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i)=0
virtual dof_id_type key() const
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 int num_sides
Geometric constants for all Tets.
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
Node * _nodelinks_data[num_nodes]
Data for links to nodes.
This is the base class from which all geometric element types are derived.
Tet4(Elem *p=nullptr)
Constructor.
static const int num_edges
virtual bool is_edge(const unsigned int i) const override
virtual unsigned int n_sub_elem() const override
The libMesh namespace provides an interface to certain functionality in the library.
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
virtual dof_id_type key() const override
The Tet is an element in 3D composed of 4 sides.
static const int nodes_per_side
Tet4 & operator=(const Tet4 &)=delete
static const int num_nodes
Geometric constants for Tet4.
virtual Order default_order() const override
static const Real _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
ElemType side_type(const unsigned int s) const override final
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Builds a TRI3 built coincident with face i.
virtual unsigned int n_nodes() const override
virtual bool is_face(const unsigned int i) const override
virtual bool has_invertible_map(Real tol) const override
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
virtual bool is_linear() const override
virtual Real volume() const override
An optimized method for computing the area of a 4-node tetrahedron.
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.
LIBMESH_ENABLE_TOPOLOGY_CACHES
The Tet4 is an element in 3D composed of 4 nodes.
static const int num_children
virtual ElemType type() const override
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
This maps the node of the edge to element node numbers.
virtual bool has_affine_map() const override
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE2 built coincident with face i.
virtual Real embedding_matrix(const unsigned int i, const unsigned int j, const unsigned int k) const override
Matrix used to create the elements children.
A Point defines a location in LIBMESH_DIM dimensional Real space.
static const int nodes_per_edge