libMesh
face_tri3.C
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 // Local includes
19 #include "libmesh/edge_edge2.h"
20 #include "libmesh/face_tri3.h"
21 #include "libmesh/enum_io_package.h"
22 #include "libmesh/enum_order.h"
23 
24 namespace libMesh
25 {
26 
27 
28 
29 // ------------------------------------------------------------
30 // Tri3 class static member initializations
31 const int Tri3::num_nodes;
32 const int Tri3::nodes_per_side;
33 
35  {
36  {0, 1}, // Side 0
37  {1, 2}, // Side 1
38  {2, 0} // Side 2
39  };
40 
41 #ifdef LIBMESH_ENABLE_AMR
42 
44  {
45  // embedding matrix for child 0
46  {
47  // 0 1 2
48  {1.0, 0.0, 0.0}, // 0
49  {0.5, 0.5, 0.0}, // 1
50  {0.5, 0.0, 0.5} // 2
51  },
52 
53  // embedding matrix for child 1
54  {
55  // 0 1 2
56  {0.5, 0.5, 0.0}, // 0
57  {0.0, 1.0, 0.0}, // 1
58  {0.0, 0.5, 0.5} // 2
59  },
60 
61  // embedding matrix for child 2
62  {
63  // 0 1 2
64  {0.5, 0.0, 0.5}, // 0
65  {0.0, 0.5, 0.5}, // 1
66  {0.0, 0.0, 1.0} // 2
67  },
68 
69  // embedding matrix for child 3
70  {
71  // 0 1 2
72  {0.5, 0.5, 0.0}, // 0
73  {0.0, 0.5, 0.5}, // 1
74  {0.5, 0.0, 0.5} // 2
75  }
76  };
77 
78 #endif
79 
80 
81 
82 // ------------------------------------------------------------
83 // Tri3 class member functions
84 
85 bool Tri3::is_vertex(const unsigned int libmesh_dbg_var(n)) const
86 {
87  libmesh_assert_not_equal_to (n, invalid_uint);
88  return true;
89 }
90 
91 bool Tri3::is_edge(const unsigned int) const
92 {
93  return false;
94 }
95 
96 bool Tri3::is_face(const unsigned int) const
97 {
98  return false;
99 }
100 
101 bool Tri3::is_node_on_side(const unsigned int n,
102  const unsigned int s) const
103 {
104  libmesh_assert_less (s, n_sides());
105  return std::find(std::begin(side_nodes_map[s]),
106  std::end(side_nodes_map[s]),
107  n) != std::end(side_nodes_map[s]);
108 }
109 
110 std::vector<unsigned>
111 Tri3::nodes_on_side(const unsigned int s) const
112 {
113  libmesh_assert_less(s, n_sides());
114  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
115 }
116 
117 std::vector<unsigned>
118 Tri3::nodes_on_edge(const unsigned int e) const
119 {
120  return nodes_on_side(e);
121 }
122 
124 {
125  return FIRST;
126 }
127 
129 {
130  return this->volume() > tol;
131 }
132 
133 std::unique_ptr<Elem> Tri3::build_side_ptr (const unsigned int i)
134 {
135  return this->simple_build_side_ptr<Edge2, Tri3>(i);
136 }
137 
138 
139 
140 void Tri3::build_side_ptr (std::unique_ptr<Elem> & side,
141  const unsigned int i)
142 {
143  this->simple_build_side_ptr<Tri3>(side, i, EDGE2);
144 }
145 
146 
147 
148 void Tri3::connectivity(const unsigned int libmesh_dbg_var(sf),
149  const IOPackage iop,
150  std::vector<dof_id_type> & conn) const
151 {
152  libmesh_assert_less (sf, this->n_sub_elem());
153  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
154 
155  switch (iop)
156  {
157  case TECPLOT:
158  {
159  conn.resize(4);
160  conn[0] = this->node_id(0)+1;
161  conn[1] = this->node_id(1)+1;
162  conn[2] = this->node_id(2)+1;
163  conn[3] = this->node_id(2)+1;
164  return;
165  }
166 
167  case VTK:
168  {
169  conn.resize(3);
170  conn[0] = this->node_id(0);
171  conn[1] = this->node_id(1);
172  conn[2] = this->node_id(2);
173  return;
174  }
175 
176  default:
177  libmesh_error_msg("Unsupported IO package " << iop);
178  }
179 }
180 
181 
182 
183 
184 
185 
187 {
188  return Elem::vertex_average();
189 }
190 
192 {
193  // 3-node triangles have the following formula for computing the area
194  return 0.5 * cross_norm(point(1) - point(0),
195  point(2) - point(0));
196 }
197 
198 
199 
200 std::pair<Real, Real> Tri3::min_and_max_angle() const
201 {
202  Point v10 ( this->point(1) - this->point(0) );
203  Point v20 ( this->point(2) - this->point(0) );
204  Point v21 ( this->point(2) - this->point(1) );
205 
206  const Real
207  len_10=v10.norm(),
208  len_20=v20.norm(),
209  len_21=v21.norm();
210 
211  const Real
212  theta0=std::acos(( v10*v20)/len_10/len_20),
213  theta1=std::acos((-v10*v21)/len_10/len_21),
214  theta2=libMesh::pi - theta0 - theta1
215  ;
216 
217  libmesh_assert_greater (theta0, 0.);
218  libmesh_assert_greater (theta1, 0.);
219  libmesh_assert_greater (theta2, 0.);
220 
221  return std::make_pair(std::min(theta0, std::min(theta1,theta2)),
222  std::max(theta0, std::max(theta1,theta2)));
223 }
224 
225 bool Tri3::contains_point (const Point & p, Real tol) const
226 {
227  // See "Barycentric Technique" section at
228  // http://www.blackpawn.com/texts/pointinpoly for details.
229 
230  // Compute vectors
231  Point v0 = this->point(1) - this->point(0);
232  Point v1 = this->point(2) - this->point(0);
233  Point v2 = p - this->point(0);
234 
235  // Compute dot products
236  Real dot00 = v0 * v0;
237  Real dot01 = v0 * v1;
238  Real dot02 = v0 * v2;
239  Real dot11 = v1 * v1;
240  Real dot12 = v1 * v2;
241 
242  // Out of plane check
243  if (std::abs(triple_product(v2, v0, v1)) / std::max(dot00, dot11) > tol)
244  return false;
245 
246  // Compute barycentric coordinates
247  Real invDenom = 1 / (dot00 * dot11 - dot01 * dot01);
248  Real u = (dot11 * dot02 - dot01 * dot12) * invDenom;
249  Real v = (dot00 * dot12 - dot01 * dot02) * invDenom;
250 
251  // Check if point is in triangle
252  return (u > -tol) && (v > -tol) && (u + v < 1 + tol);
253 }
254 
257 {
258  return Elem::loose_bounding_box();
259 }
260 
261 
262 void Tri3::permute(unsigned int perm_num)
263 {
264  libmesh_assert_less (perm_num, 3);
265 
266  for (unsigned int i = 0; i != perm_num; ++i)
267  {
268  swap3nodes(0,1,2);
269  swap3neighbors(0,1,2);
270  }
271 }
272 
273 
274 void Tri3::flip(BoundaryInfo * boundary_info)
275 {
276  libmesh_assert(boundary_info);
277 
278  swap2nodes(0,1);
279  swap2neighbors(1,2);
280  swap2boundarysides(1,2,boundary_info);
281  swap2boundaryedges(1,2,boundary_info);
282 }
283 
284 
285 ElemType
286 Tri3::side_type (const unsigned int libmesh_dbg_var(s)) const
287 {
288  libmesh_assert_less (s, 3);
289  return EDGE2;
290 }
291 
292 } // namespace libMesh
ElemType
Defines an enum for geometric element types.
void swap2boundaryedges(unsigned short e1, unsigned short e2, BoundaryInfo *boundary_info) const
Swaps two edges in boundary_info, if it is non-null.
Definition: elem.C:3550
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
auto norm() const -> decltype(std::norm(T()))
Definition: type_vector.h:907
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
static const int num_nodes
Geometric constants for Tri3.
Definition: face_tri3.h:162
virtual unsigned int n_sub_elem() const override
Definition: face_tri3.h:86
static const int num_sides
Geometric constants for all Tris.
Definition: face_tri.h:86
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
Definition: face_tri3.h:169
virtual Point true_centroid() const override
The centroid of a 3-node triangle is simply given by the average of its vertex positions.
Definition: face_tri3.C:186
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Calls cross_norm_sq() and takes the square root of the result.
Definition: type_vector.h:1096
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: face_tri3.C:101
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
virtual BoundingBox loose_bounding_box() const
Definition: elem.C:3465
std::pair< Real, Real > min_and_max_angle() const
Definition: face_tri3.C:200
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: face_tri3.C:111
void swap2boundarysides(unsigned short s1, unsigned short s2, BoundaryInfo *boundary_info) const
Swaps two sides in boundary_info, if it is non-null.
Definition: elem.C:3534
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Definition: face_tri3.C:133
virtual BoundingBox loose_bounding_box() const override
Builds a bounding box out of the nodal positions.
Definition: face_tri3.C:256
virtual bool has_invertible_map(Real tol) const override
Definition: face_tri3.C:128
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
Definition: face_tri.h:87
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: face_tri3.C:148
ElemType side_type(const unsigned int s) const override final
Definition: face_tri3.C:286
virtual bool is_edge(const unsigned int i) const override
Definition: face_tri3.C:91
void swap3neighbors(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three neighbor_ptrs, "rotating" them.
Definition: elem.h:2133
virtual unsigned int n_sides() const override final
Definition: face_tri.h:98
void swap3nodes(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three node_ptrs, "rotating" them.
Definition: elem.h:2124
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1029
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
virtual bool contains_point(const Point &p, Real tol) const override
Specialization for tri3 elements.
Definition: face_tri3.C:225
static const int nodes_per_side
Definition: face_tri3.h:163
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:2102
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
virtual bool is_vertex(const unsigned int i) const override
Definition: face_tri3.C:85
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...
Definition: face_tri3.C:262
virtual Order default_order() const override
Definition: face_tri3.C:123
static const Real _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
Definition: face_tri3.h:229
virtual bool is_face(const unsigned int i) const override
Definition: face_tri3.C:96
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
Definition: face_tri3.C:274
virtual Real volume() const override
An optimized method for computing the area of a 3-node triangle.
Definition: face_tri3.C:191
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:2475
const Point & point(const unsigned int i) const
Definition: elem.h:2453
Point vertex_average() const
Definition: elem.C:688
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
Definition: face_tri3.C:118
const Real pi
.
Definition: libmesh.h:299