face_tri3.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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_TRI3_H
21 #define LIBMESH_FACE_TRI3_H
22 
23 
24 // Local includes
25 #include "libmesh/libmesh_common.h"
26 #include "libmesh/face_tri.h"
27 
28 // C++ includes
29 #include <cstddef>
30 
31 namespace libMesh
32 {
33 
56 class Tri3 : public Tri
57 {
58 public:
59 
63  explicit
64  Tri3 (Elem * p=nullptr) :
65  Tri(Tri3::n_nodes(), p, _nodelinks_data) {}
66 
67  Tri3 (Tri3 &&) = delete;
68  Tri3 (const Tri3 &) = delete;
69  Tri3 & operator= (const Tri3 &) = delete;
70  Tri3 & operator= (Tri3 &&) = delete;
71  virtual ~Tri3() = default;
72 
76  virtual ElemType type () const override { return TRI3; }
77 
81  virtual unsigned int n_sub_elem() const override { return 1; }
82 
86  virtual bool is_vertex(const unsigned int i) const override;
87 
91  virtual bool is_edge(const unsigned int i) const override;
92 
96  virtual bool is_face(const unsigned int i) const override;
97 
102  virtual bool is_node_on_side(const unsigned int n,
103  const unsigned int s) const override;
104 
105  virtual std::vector<unsigned int> nodes_on_side(const unsigned int s) const override;
106 
111  virtual bool is_node_on_edge(const unsigned int n,
112  const unsigned int e) const override
113  { return this->is_node_on_side(n,e); }
114 
119  virtual bool has_affine_map () const override { return true; }
120 
125  virtual bool is_linear () const override { return true; }
126 
130  virtual Order default_order() const override;
131 
132  virtual std::unique_ptr<Elem> build_side_ptr (const unsigned int i,
133  bool proxy) override;
134 
138  virtual void build_side_ptr (std::unique_ptr<Elem> & elem,
139  const unsigned int i) override;
140 
141  virtual void connectivity(const unsigned int sf,
142  const IOPackage iop,
143  std::vector<dof_id_type> & conn) const override;
144 
148  static const int num_nodes = 3;
149  static const int num_sides = 3;
150  static const int num_children = 4;
151  static const int nodes_per_side = 2;
152 
157  static const unsigned int side_nodes_map[num_sides][nodes_per_side];
158 
162  virtual Real volume () const override;
163 
169  std::pair<Real, Real> min_and_max_angle() const;
170 
175  virtual bool contains_point (const Point & p, Real tol) const override;
176 
180  virtual BoundingBox loose_bounding_box () const override;
181 
182 protected:
183 
188 
189 
190 
191 #ifdef LIBMESH_ENABLE_AMR
192 
196  virtual float embedding_matrix (const unsigned int i,
197  const unsigned int j,
198  const unsigned int k) const override
199  { return _embedding_matrix[i][j][k]; }
200 
206 
208 
209 #endif // LIBMESH_ENABLE_AMR
210 
211 };
212 
213 } // namespace libMesh
214 
215 #endif // LIBMESH_FACE_TRI3_H
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: face_tri3.h:111
A 2D triangular element with 3 nodes.
Definition: face_tri3.h:56
Node * _nodelinks_data[num_nodes]
Definition: face_tri3.h:187
A geometric point in (x,y,z) space associated with a DOF.
Definition: node.h:52
static const int num_nodes
Definition: face_tri3.h:148
virtual unsigned int n_sub_elem() const override
Definition: face_tri3.h:81
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
Definition: face_tri3.h:157
static const int num_sides
Definition: face_tri3.h:149
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: face_tri3.C:104
The base class for all geometric element types.
Definition: elem.h:100
The base class for all triangular element types.
Definition: face_tri.h:47
std::pair< Real, Real > min_and_max_angle() const
Definition: face_tri3.C:203
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: face_tri3.C:114
virtual BoundingBox loose_bounding_box() const override
Definition: face_tri3.C:259
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: face_tri3.C:156
virtual bool is_edge(const unsigned int i) const override
Definition: face_tri3.C:94
virtual ElemType type() const override
Definition: face_tri3.h:76
virtual bool contains_point(const Point &p, Real tol) const override
Definition: face_tri3.C:228
virtual bool is_linear() const override
Definition: face_tri3.h:125
LIBMESH_ENABLE_TOPOLOGY_CACHES
Definition: face_tri3.h:207
static const int nodes_per_side
Definition: face_tri3.h:151
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy) override
Definition: face_tri3.C:125
virtual unsigned int n_nodes() const override
Definition: face_tri.h:87
virtual float embedding_matrix(const unsigned int i, const unsigned int j, const unsigned int k) const override
Definition: face_tri3.h:196
virtual bool is_vertex(const unsigned int i) const override
Definition: face_tri3.C:89
Tri3 & operator=(const Tri3 &)=delete
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Order default_order() const override
Definition: face_tri3.C:120
virtual bool has_affine_map() const override
Definition: face_tri3.h:119
static const int num_children
Definition: face_tri3.h:150
virtual bool is_face(const unsigned int i) const override
Definition: face_tri3.C:99
virtual ~Tri3()=default
virtual Real volume() const override
Definition: face_tri3.C:194
A geometric point in (x,y,z) space.
Definition: point.h:38
static const float _embedding_matrix[num_children][num_nodes][num_nodes]
Definition: face_tri3.h:205
Tri3(Elem *p=nullptr)
Definition: face_tri3.h:64