libMesh
face_tri.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_FACE_TRI_H
21 #define LIBMESH_FACE_TRI_H
22 
23 // Local includes
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/face.h"
26 
27 namespace libMesh
28 {
29 
47 class Tri : public Face
48 {
49 public:
50 
55  Tri (const unsigned int nn,
56  Elem * p,
57  Node ** nodelinkdata) :
58  Face(nn, num_sides, p, _elemlinks_data, nodelinkdata)
59  {
60  // Make sure the interior parent isn't undefined
61  if (LIBMESH_DIM > 2)
62  this->set_interior_parent(nullptr);
63  }
64 
65  Tri (Tri &&) = delete;
66  Tri (const Tri &) = delete;
67  Tri & operator= (const Tri &) = delete;
68  Tri & operator= (Tri &&) = delete;
69  virtual ~Tri() = default;
70 
75  virtual Point master_point (const unsigned int i) const override final
76  {
77  libmesh_assert_less(i, this->n_nodes());
78  return Point(_master_points[i][0],
79  _master_points[i][1],
80  _master_points[i][2]);
81  }
82 
86  static const int num_sides = 3;
87  static const int num_children = 4;
88 
93  virtual unsigned int n_nodes() const override { return 3; }
94 
98  virtual unsigned int n_sides() const override final { return 3; }
99 
103  virtual unsigned int n_vertices() const override final { return 3; }
104 
108  virtual unsigned int n_edges() const override final { return 3; }
109 
113  virtual unsigned int n_children() const override final { return 4; }
114 
119  virtual bool is_child_on_side(const unsigned int c,
120  const unsigned int s) const override final;
121 
125  using Elem::key;
126 
132  virtual dof_id_type key () const override final;
133 
138  virtual dof_id_type key (const unsigned int s) const override;
139 
146  virtual dof_id_type low_order_key (const unsigned int s) const override;
147 
151  virtual unsigned int local_side_node(unsigned int side,
152  unsigned int side_node) const override;
153 
159  virtual unsigned int local_edge_node(unsigned int edge,
160  unsigned int edge_node) const override;
161 
165  virtual std::unique_ptr<Elem> side_ptr (const unsigned int i) override final;
166 
170  virtual void side_ptr (std::unique_ptr<Elem> & elem,
171  const unsigned int i) override final;
172 
177  virtual Real quality (const ElemQuality q) const override;
178 
184  virtual std::pair<Real, Real> qual_bounds (const ElemQuality q) const override;
185 
189  virtual unsigned int n_permutations() const override final { return 3; }
190 
191  virtual bool is_flipped() const override final;
192 
193  virtual std::vector<unsigned int> edges_adjacent_to_node(const unsigned int n) const override;
194 
195  virtual bool on_reference_element(const Point & p,
196  const Real eps = TOLERANCE) const override final;
197 
198 protected:
199 
203  Elem * _elemlinks_data[4+(LIBMESH_DIM>2)];
204 
208  static const Real _master_points[6][3];
209 
217  static const unsigned int adjacent_sides_map[/*num_vertices*/3][/*n_adjacent_sides*/2];
218 };
219 
220 } // namespace libMesh
221 
222 #endif // LIBMESH_FACE_TRI_H
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
Definition: face_tri.C:107
static const Real _master_points[6][3]
Master element node locations.
Definition: face_tri.h:208
Tri(const unsigned int nn, Elem *p, Node **nodelinkdata)
Default triangular element, takes number of nodes and parent.
Definition: face_tri.h:55
A Node is like a Point, but with more information.
Definition: node.h:52
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
Calls local_side_node(edge, edge_node).
Definition: face_tri.C:90
The Face is an abstract element type that lives in two dimensions.
Definition: face.h:37
static const int num_sides
Geometric constants for all Tris.
Definition: face_tri.h:86
virtual dof_id_type key() const
Definition: elem.C:753
static constexpr Real TOLERANCE
virtual bool is_flipped() const override final
Definition: face_tri.C:139
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
The Tri is an element in 2D composed of 3 sides.
Definition: face_tri.h:47
virtual unsigned int n_vertices() const override final
Definition: face_tri.h:103
virtual unsigned int n_permutations() const override final
Three sides, one orientation each.
Definition: face_tri.h:189
virtual ~Tri()=default
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
Definition: face_tri.h:87
void set_interior_parent(Elem *p)
Sets the pointer to the element&#39;s interior_parent.
Definition: elem.C:1248
virtual unsigned int n_sides() const override final
Definition: face_tri.h:98
Tri & operator=(const Tri &)=delete
virtual std::vector< unsigned int > edges_adjacent_to_node(const unsigned int n) const override
Definition: face_tri.C:155
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override final
Definition: face_tri.C:129
virtual dof_id_type key() const override final
Definition: face_tri.C:98
virtual dof_id_type low_order_key(const unsigned int s) const override
Definition: face_tri.C:69
ElemQuality
Defines an enum for element quality metrics.
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Definition: face_tri.C:79
virtual Point master_point(const unsigned int i) const override final
Definition: face_tri.h:75
virtual unsigned int n_nodes() const override
Definition: face_tri.h:93
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
Definition: face_tri.C:375
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
Definition: face_tri.C:315
virtual unsigned int n_children() const override final
Definition: face_tri.h:113
Elem * _elemlinks_data[4+(LIBMESH_DIM >2)]
Data for links to parent/neighbor/interior_parent elements.
Definition: face_tri.h:203
virtual Real quality(const ElemQuality q) const override
Definition: face_tri.C:172
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
virtual unsigned int n_edges() const override final
Definition: face_tri.h:108
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.
Definition: face_tri.h:217
uint8_t dof_id_type
Definition: id_types.h:67