20 #ifndef LIBMESH_FACE_TRI3_H 21 #define LIBMESH_FACE_TRI3_H 25 #include "libmesh/libmesh_common.h" 26 #include "libmesh/face_tri.h" 76 virtual ~Tri3() =
default;
86 virtual unsigned int n_sub_elem()
const override {
return 1; }
91 virtual bool is_vertex(
const unsigned int i)
const override;
96 virtual bool is_edge(
const unsigned int i)
const override;
101 virtual bool is_face(
const unsigned int i)
const override;
108 const unsigned int s)
const override;
110 virtual std::vector<unsigned int>
nodes_on_side(
const unsigned int s)
const override;
112 virtual std::vector<unsigned int>
nodes_on_edge(
const unsigned int e)
const override;
119 const unsigned int e)
const override 137 virtual bool is_linear ()
const override {
return true; }
144 virtual std::unique_ptr<Elem>
build_side_ptr (
const unsigned int i)
override;
150 const unsigned int i)
override;
157 std::vector<dof_id_type> & conn)
const override;
200 virtual void permute(
unsigned int perm_num)
override final;
215 #ifdef LIBMESH_ENABLE_AMR 221 const unsigned int j,
222 const unsigned int k)
const override 233 #endif // LIBMESH_ENABLE_AMR 239 #endif // LIBMESH_FACE_TRI3_H virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
The Tri3 is an element in 2D composed of 3 nodes.
ElemType
Defines an enum for geometric element types.
Node * _nodelinks_data[num_nodes]
Data for links to nodes.
Order
defines an enum for polynomial orders.
A Node is like a Point, but with more information.
static const int num_nodes
Geometric constants for Tri3.
virtual unsigned int n_sub_elem() const override
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.
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i)=0
virtual Point true_centroid() const override
The centroid of a 3-node triangle is simply given by the average of its vertex positions.
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
This is the base class from which all geometric element types are derived.
The Tri is an element in 2D composed of 3 sides.
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) override
virtual BoundingBox loose_bounding_box() const override
Builds a bounding box out of the nodal positions.
virtual bool has_invertible_map(Real tol) const override
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
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.
ElemType side_type(const unsigned int s) const override final
virtual bool is_edge(const unsigned int i) const override
virtual ElemType type() const override
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
virtual bool contains_point(const Point &p, Real tol) const override
Specialization for tri3 elements.
virtual bool is_linear() const override
LIBMESH_ENABLE_TOPOLOGY_CACHES
static const int nodes_per_side
Defines a Cartesian bounding box by the two corner extremum.
virtual bool is_vertex(const unsigned int i) const override
Tri3 & operator=(const Tri3 &)=delete
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void permute(unsigned int perm_num) override final
Permutes the element (by swapping node and neighbor pointers) according to the specified index...
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.
virtual bool has_affine_map() const override
virtual bool is_face(const unsigned int i) const override
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
virtual Real volume() const override
An optimized method for computing the area of a 3-node triangle.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Tri3(Elem *p=nullptr)
Constructor.
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override