libMesh
face_tri7.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_edge3.h"
20 #include "libmesh/face_tri7.h"
21 #include "libmesh/enum_io_package.h"
22 #include "libmesh/enum_order.h"
23 
24 #ifdef LIBMESH_ENABLE_AMR
25 namespace {
26  constexpr libMesh::Real r18 = 18;
27 }
28 #endif
29 
30 namespace libMesh
31 {
32 
33 
34 
35 
36 // ------------------------------------------------------------
37 // Tri7 class static member initializations
38 const int Tri7::num_nodes;
39 const int Tri7::nodes_per_side;
40 
42  {
43  {0, 1, 3}, // Side 0
44  {1, 2, 4}, // Side 1
45  {2, 0, 5} // Side 2
46  };
47 
48 
49 #ifdef LIBMESH_ENABLE_AMR
50 
52  {
53  // embedding matrix for child 0
54  {
55  // 0 1 2 3 4 5 6
56  { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0
57  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 1
58  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 2
59  { .375, -.125, 0.0, .75, 0.0, 0.0, 0.0}, // 3
60  {.09375,-.03125,-.03125, .125, -.125, .125,.84375}, // 4
61  { .375, 0.0, -.125, 0.0, 0.0, .75, 0.0}, // 5
62  { 5/r18,-1/r18,-1/r18,4/r18,-2/r18,4/r18, 0.5} // 6
63  },
64 
65  // embedding matrix for child 1
66  {
67  // 0 1 2 3 4 5 6
68  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 0
69  { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 1
70  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 2
71  { -.125, .375, 0.0, .75, 0.0, 0.0, 0.0}, // 3
72  { 0.0, .375, -.125, 0.0, .75, 0.0, 0.0}, // 4
73  {-.03125,.09375,-.03125, .125, .125, -.125,.84375}, // 5
74  {-1/r18, 5/r18,-1/r18,4/r18,4/r18,-2/r18, 0.5} // 6
75  },
76 
77  // embedding matrix for child 2
78  {
79  // 0 1 2 3 4 5 6
80  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 0
81  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 1
82  { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 2
83  {-.03125,-.03125,.09375, -.125, .125, .125,.84375}, // 3
84  { 0.0, -.125, .375, 0.0, .75, 0.0, 0.0}, // 4
85  { -.125, 0.0, .375, 0.0, 0.0, .75, 0.0}, // 5
86  {-1/r18,-1/r18, 5/r18,-2/r18,4/r18,4/r18, 0.5} // 6
87  },
88 
89  // embedding matrix for child 3
90  {
91  // 0 1 2 3 4 5 6
92  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 0
93  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 1
94  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 2
95  {-.03125,.09375,-.03125, .125, .125,-.125,.84375}, // 3
96  {-.03125,-.03125,.09375,-.125, .125, .125,.84375}, // 4
97  {.09375,-.03125,-.03125, .125,-.125, .125,.84375}, // 5
98  { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0} // 6
99  }
100  };
101 
102 const std::vector<std::pair<unsigned char, unsigned char>>
104  {
105  // Child 0
106  { {},{{0,1}},{{0,2}},{{0,3}},{{3,5}},{{0,5}},{{0,6}} },
107  // Child 1
108  { {{0,1}}, {},{{1,2}},{{1,3}},{{1,4}},{{3,4}},{{1,6}} },
109  // Child 2
110  { {{0,2}},{{1,2}}, {},{{4,5}},{{2,4}},{{2,5}},{{2,6}} },
111  // Child 3
112  { {{0,1}},{{1,2}},{{0,2}},{{3,4}},{{4,5}},{{3,5}},{{0,4},{1,5},{2,3}} }
113  };
114 #endif
115 
116 
117 
118 // ------------------------------------------------------------
119 // Tri7 class member functions
120 
121 bool Tri7::is_vertex(const unsigned int i) const
122 {
123  if (i < 3)
124  return true;
125  return false;
126 }
127 
128 bool Tri7::is_edge(const unsigned int i) const
129 {
130  if (i < 3 || i == 6)
131  return false;
132  return true;
133 }
134 
135 bool Tri7::is_face(const unsigned int i) const
136 {
137  if (i > 5)
138  return true;
139  return false;
140 }
141 
142 bool Tri7::is_node_on_side(const unsigned int n,
143  const unsigned int s) const
144 {
145  libmesh_assert_less (s, n_sides());
146  return std::find(std::begin(side_nodes_map[s]),
147  std::end(side_nodes_map[s]),
148  n) != std::end(side_nodes_map[s]);
149 }
150 
151 std::vector<unsigned>
152 Tri7::nodes_on_side(const unsigned int s) const
153 {
154  libmesh_assert_less(s, n_sides());
155  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
156 }
157 
158 std::vector<unsigned>
159 Tri7::nodes_on_edge(const unsigned int e) const
160 {
161  return nodes_on_side(e);
162 }
163 
165 {
166  // Make sure edges are straight
167  Point v = this->point(2) - this->point(1);
168  if (!v.relative_fuzzy_equals
169  ((this->point(4) - this->point(1))*2))
170  return false;
171  v = this->point(1) - this->point(0);
172  if (!v.relative_fuzzy_equals
173  ((this->point(3) - this->point(0))*2))
174  return false;
175  Point v20 = this->point(2) - this->point(0);
176  if (!v20.relative_fuzzy_equals
177  ((this->point(5) - this->point(0))*2))
178  return false;
179 
180  // Make sure center node is centered
181  v += v20;
182  if (!v.relative_fuzzy_equals
183  ((this->point(6) - this->point(0))*3))
184  return false;
185 
186  return true;
187 }
188 
189 
190 
192 {
193  return THIRD;
194 }
195 
196 
197 
199 {
200  return SECOND;
201 }
202 
203 
204 
205 dof_id_type Tri7::key (const unsigned int s) const
206 {
207  libmesh_assert_less (s, this->n_sides());
208 
209  switch (s)
210  {
211  case 0:
212 
213  return
214  this->compute_key (this->node_id(3));
215 
216  case 1:
217 
218  return
219  this->compute_key (this->node_id(4));
220 
221  case 2:
222 
223  return
224  this->compute_key (this->node_id(5));
225 
226  default:
227  libmesh_error_msg("Invalid side s = " << s);
228  }
229 }
230 
231 
232 
233 unsigned int Tri7::local_side_node(unsigned int side,
234  unsigned int side_node) const
235 {
236  libmesh_assert_less (side, this->n_sides());
237  libmesh_assert_less (side_node, Tri7::nodes_per_side);
238 
239  return Tri7::side_nodes_map[side][side_node];
240 }
241 
242 
243 
244 std::unique_ptr<Elem> Tri7::build_side_ptr (const unsigned int i)
245 {
246  return this->simple_build_side_ptr<Edge3, Tri7>(i);
247 }
248 
249 
250 
251 void Tri7::build_side_ptr (std::unique_ptr<Elem> & side,
252  const unsigned int i)
253 {
254  this->simple_build_side_ptr<Tri7>(side, i, EDGE3);
255 }
256 
257 
258 
259 void Tri7::connectivity(const unsigned int sf,
260  const IOPackage iop,
261  std::vector<dof_id_type> & conn) const
262 {
263  libmesh_assert_less (sf, this->n_sub_elem());
264  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
265 
266  switch (iop)
267  {
268  case TECPLOT:
269  {
270  conn.resize(4);
271  switch(sf)
272  {
273  case 0:
274  // linear sub-triangle 0
275  conn[0] = this->node_id(0)+1;
276  conn[1] = this->node_id(3)+1;
277  conn[2] = this->node_id(5)+1;
278  conn[3] = this->node_id(5)+1;
279 
280  return;
281 
282  case 1:
283  // linear sub-triangle 1
284  conn[0] = this->node_id(3)+1;
285  conn[1] = this->node_id(1)+1;
286  conn[2] = this->node_id(4)+1;
287  conn[3] = this->node_id(4)+1;
288 
289  return;
290 
291  case 2:
292  // linear sub-triangle 2
293  conn[0] = this->node_id(5)+1;
294  conn[1] = this->node_id(4)+1;
295  conn[2] = this->node_id(2)+1;
296  conn[3] = this->node_id(2)+1;
297 
298  return;
299 
300  case 3:
301  // linear sub-triangle 3
302  conn[0] = this->node_id(3)+1;
303  conn[1] = this->node_id(4)+1;
304  conn[2] = this->node_id(5)+1;
305  conn[3] = this->node_id(5)+1;
306 
307  return;
308 
309  default:
310  libmesh_error_msg("Invalid sf = " << sf);
311  }
312  }
313 
314  case VTK:
315  {
316  // VTK has a vtkBiQuadraticTriangle class whose connectivity matches libMesh's
317  conn.resize(Tri7::num_nodes);
318  for (auto i : index_range(conn))
319  conn[i] = this->node_id(i);
320  return;
321  }
322 
323  default:
324  libmesh_error_msg("Unsupported IO package " << iop);
325  }
326 }
327 
328 
329 
331 {
332  // This might have curved edges, or might be a curved surface in
333  // 3-space, in which case the full bounding box can be larger than
334  // the bounding box of just the nodes.
335  //
336  //
337  // FIXME - I haven't yet proven the formula below to be correct for
338  // quadratics in 2D - RHS
339  //
340  // FIXME - This doesn't take into account curvature caused by the
341  // center node in 3D - RHS
342  Point pmin, pmax;
343 
344  for (unsigned d=0; d<LIBMESH_DIM; ++d)
345  {
346  Real center = this->point(0)(d);
347  for (unsigned int p=1; p != 6; ++p)
348  center += this->point(p)(d);
349  center /= 6;
350 
351  Real hd = std::abs(center - this->point(0)(d));
352  for (unsigned int p=1; p != 6; ++p)
353  hd = std::max(hd, std::abs(center - this->point(p)(d)));
354 
355  pmin(d) = center - hd;
356  pmax(d) = center + hd;
357  }
358 
359  return BoundingBox(pmin, pmax);
360 }
361 
362 
363 
364 unsigned int Tri7::n_second_order_adjacent_vertices (const unsigned int n) const
365 {
366  switch (n)
367  {
368  case 3:
369  case 4:
370  case 5:
371  return 2;
372 
373  case 6:
374  return 3;
375 
376  default:
377  libmesh_error_msg("Invalid n = " << n);
378  }
379 }
380 
381 
382 
383 unsigned short int Tri7::second_order_adjacent_vertex (const unsigned int n,
384  const unsigned int v) const
385 {
386  libmesh_assert_greater_equal (n, this->n_vertices());
387  libmesh_assert_less (n, this->n_nodes());
388 
389  switch (n)
390  {
391  case 6:
392  {
393  libmesh_assert_less (v, 3);
394  return static_cast<unsigned short int>(v);
395  }
396 
397  default:
398  {
399  libmesh_assert_less (v, 2);
400  return _second_order_adjacent_vertices[n-this->n_vertices()][v];
401  }
402  }
403 }
404 
405 
406 
407 const unsigned short int Tri7::_second_order_adjacent_vertices[Tri7::num_sides][2] =
408  {
409  {0, 1}, // vertices adjacent to node 3
410  {1, 2}, // vertices adjacent to node 4
411  {0, 2} // vertices adjacent to node 5
412  };
413 
414 
415 
416 std::pair<unsigned short int, unsigned short int>
417 Tri7::second_order_child_vertex (const unsigned int n) const
418 {
419  libmesh_assert_greater_equal (n, this->n_vertices());
420  libmesh_assert_less (n, this->n_nodes());
421  return std::pair<unsigned short int, unsigned short int>
424 }
425 
426 
427 
429  {
430  99,99,99, // Vertices
431  0,1,0, // Edges
432  3 // Interior
433  };
434 
435 
436 
438  {
439  99,99,99, // Vertices
440  1,2,2, // Edges
441  6 // Interior
442  };
443 
444 
445 void Tri7::permute(unsigned int perm_num)
446 {
447  libmesh_assert_less (perm_num, 3);
448 
449  for (unsigned int i = 0; i != perm_num; ++i)
450  {
451  swap3nodes(0,1,2);
452  swap3nodes(3,4,5);
453  swap3neighbors(0,1,2);
454  }
455 }
456 
457 
458 void Tri7::flip(BoundaryInfo * boundary_info)
459 {
460  libmesh_assert(boundary_info);
461 
462  swap2nodes(0,1);
463  swap2nodes(4,5);
464  swap2neighbors(1,2);
465  swap2boundarysides(1,2,boundary_info);
466  swap2boundaryedges(1,2,boundary_info);
467 }
468 
469 
470 unsigned int Tri7::center_node_on_side(const unsigned short side) const
471 {
472  libmesh_assert_less (side, Tri7::num_sides);
473  return side + 3;
474 }
475 
476 
477 ElemType
478 Tri7::side_type (const unsigned int libmesh_dbg_var(s)) const
479 {
480  libmesh_assert_less (s, 3);
481  return EDGE3;
482 }
483 
484 
485 } // namespace libMesh
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: face_tri7.C:152
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
Definition: face_tri7.C:417
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
virtual bool has_affine_map() const override
Definition: face_tri7.C:164
static const unsigned short int _second_order_adjacent_vertices[num_sides][2]
Matrix that tells which vertices define the location of mid-side (or second-order) nodes...
Definition: face_tri7.h:277
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: face_tri7.C:142
virtual unsigned int n_second_order_adjacent_vertices(const unsigned int) const override
Definition: face_tri7.C:364
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_tri7.C:445
virtual Order default_order() const override
Definition: face_tri7.C:191
static const int num_sides
Geometric constants for all Tris.
Definition: face_tri.h:86
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Definition: face_tri7.C:244
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
virtual unsigned int n_nodes() const override
Definition: face_tri7.h:86
virtual bool is_face(const unsigned int i) const override
Definition: face_tri7.C:135
virtual unsigned int n_vertices() const override final
Definition: face_tri.h:103
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
static const unsigned short int _second_order_vertex_child_number[num_nodes]
Vector that names a child sharing each second order node.
Definition: face_tri7.h:282
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
Definition: face_tri.h:87
void swap3neighbors(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three neighbor_ptrs, "rotating" them.
Definition: elem.h:2133
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
Definition: face_tri7.C:159
static const unsigned short int _second_order_vertex_child_index[num_nodes]
Vector that names the child vertex index for each second order node.
Definition: face_tri7.h:287
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
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
virtual dof_id_type key() const override final
Definition: face_tri.C:98
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_tri7.h:258
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
virtual Order default_side_order() const override
Definition: face_tri7.C:198
static const int num_nodes
Geometric constants for Tri7.
Definition: face_tri7.h:204
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:2102
virtual bool is_vertex(const unsigned int i) const override
Definition: face_tri7.C:121
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: face_tri7.C:259
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
ElemType side_type(const unsigned int s) const override final
Definition: face_tri7.C:478
static const std::vector< std::pair< unsigned char, unsigned char > > _parent_bracketing_nodes[num_children][num_nodes]
Pairs of nodes that bracket child nodes when doing mesh refinement.
Definition: face_tri7.h:265
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned int n_sub_elem() const override
Definition: face_tri7.h:91
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_tri7.C:458
unsigned int center_node_on_side(const unsigned short side) const override final
Definition: face_tri7.C:470
virtual BoundingBox loose_bounding_box() const override
Definition: face_tri7.C:330
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:3294
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Definition: face_tri7.C:233
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
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:981
static const int nodes_per_side
Definition: face_tri7.h:205
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:117
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const override
Definition: face_tri7.C:383
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_tri7.h:211
uint8_t dof_id_type
Definition: id_types.h:67
virtual bool is_edge(const unsigned int i) const override
Definition: face_tri7.C:128