libMesh
cell_tet14.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 
19 // Local includes
20 #include "libmesh/cell_tet14.h"
21 #include "libmesh/edge_edge3.h"
22 #include "libmesh/face_tri7.h"
23 #include "libmesh/enum_io_package.h"
24 #include "libmesh/enum_order.h"
25 
26 #ifdef LIBMESH_ENABLE_AMR
27 namespace {
28  constexpr libMesh::Real r18 = 18;
29  constexpr libMesh::Real r64 = 64;
30  constexpr libMesh::Real r72 = 72;
31 }
32 #endif
33 
34 namespace libMesh
35 {
36 
37 
38 
39 // ------------------------------------------------------------
40 // Tet14 class static member initializations
41 const int Tet14::num_nodes;
42 const int Tet14::nodes_per_side;
43 const int Tet14::nodes_per_edge;
44 
46  {
47  {0, 2, 1, 6, 5, 4, 10}, // Side 0
48  {0, 1, 3, 4, 8, 7, 11}, // Side 1
49  {1, 2, 3, 5, 9, 8, 12}, // Side 2
50  {2, 0, 3, 6, 7, 9, 13} // Side 3
51  };
52 
54  {
55  {0, 1, 4}, // Edge 0
56  {1, 2, 5}, // Edge 1
57  {0, 2, 6}, // Edge 2
58  {0, 3, 7}, // Edge 3
59  {1, 3, 8}, // Edge 4
60  {2, 3, 9} // Edge 5
61  };
62 
63 // ------------------------------------------------------------
64 // Tet14 class member functions
65 
66 bool Tet14::is_vertex(const unsigned int i) const
67 {
68  if (i < 4)
69  return true;
70  return false;
71 }
72 
73 bool Tet14::is_edge(const unsigned int i) const
74 {
75  if (i < 4 || i > 9)
76  return false;
77  return true;
78 }
79 
80 bool Tet14::is_face(const unsigned int i) const
81 {
82  if (i > 9)
83  return true;
84  return false;
85 }
86 
87 bool Tet14::is_node_on_side(const unsigned int n,
88  const unsigned int s) const
89 {
90  libmesh_assert_less (s, n_sides());
91  return std::find(std::begin(side_nodes_map[s]),
92  std::end(side_nodes_map[s]),
93  n) != std::end(side_nodes_map[s]);
94 }
95 
96 std::vector<unsigned>
97 Tet14::nodes_on_side(const unsigned int s) const
98 {
99  libmesh_assert_less(s, n_sides());
100  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
101 }
102 
103 std::vector<unsigned>
104 Tet14::nodes_on_edge(const unsigned int e) const
105 {
106  libmesh_assert_less(e, n_edges());
107  return {std::begin(edge_nodes_map[e]), std::end(edge_nodes_map[e])};
108 }
109 
110 bool Tet14::is_node_on_edge(const unsigned int n,
111  const unsigned int e) const
112 {
113  libmesh_assert_less (e, n_edges());
114  return std::find(std::begin(edge_nodes_map[e]),
115  std::end(edge_nodes_map[e]),
116  n) != std::end(edge_nodes_map[e]);
117 }
118 
119 
120 #ifdef LIBMESH_ENABLE_AMR
121 
122 // This function only works if LIBMESH_ENABLE_AMR...
123 bool Tet14::is_child_on_side(const unsigned int c,
124  const unsigned int s) const
125 {
126  // Table of local IDs for the midege nodes on the side opposite a given node.
127  // See the ASCII art in the header file for this class to confirm this.
128  const unsigned int midedge_nodes_opposite[4][3] =
129  {
130  {5,8,9}, // midedge nodes opposite node 0
131  {6,7,9}, // midedge nodes opposite node 1
132  {4,7,8}, // midedge nodes opposite node 2
133  {4,5,6} // midedge nodes opposite node 3
134  };
135 
136  // Call the base class helper function
137  return Tet::is_child_on_side_helper(c, s, midedge_nodes_opposite);
138 }
139 
140 #else
141 
142 bool Tet14::is_child_on_side(const unsigned int /*c*/,
143  const unsigned int /*s*/) const
144 {
145  libmesh_not_implemented();
146  return false;
147 }
148 
149 #endif //LIBMESH_ENABLE_AMR
150 
151 
152 
154 {
155  // Make sure edges are straight
156  Point v10 = this->point(1) - this->point(0);
157  if (!v10.relative_fuzzy_equals
158  ((this->point(4) - this->point(0))*2))
159  return false;
160  Point v21 = this->point(2) - this->point(1);
161  if (!v21.relative_fuzzy_equals
162  ((this->point(5) - this->point(1))*2))
163  return false;
164  Point v20 = this->point(2) - this->point(0);
165  if (!v20.relative_fuzzy_equals
166  ((this->point(6) - this->point(0))*2))
167  return false;
168  Point v30 = this->point(3) - this->point(0);
169  if (!v30.relative_fuzzy_equals
170  ((this->point(7) - this->point(0))*2))
171  return false;
172  Point v31 = this->point(3) - this->point(1);
173  if (!v31.relative_fuzzy_equals
174  ((this->point(8) - this->point(1))*2))
175  return false;
176  Point v32 = this->point(3) - this->point(2);
177  if (!v32.relative_fuzzy_equals
178  ((this->point(9) - this->point(2))*2))
179  return false;
180 
181  // Make sure midface nodes are midface
182  if (!(v20+v10).relative_fuzzy_equals
183  ((this->point(10) - this->point(0))*3))
184  return false;
185  if (!(v30+v10).relative_fuzzy_equals
186  ((this->point(11) - this->point(0))*3))
187  return false;
188  if (!(v31+v21).relative_fuzzy_equals
189  ((this->point(12) - this->point(1))*3))
190  return false;
191  if (!(v30+v20).relative_fuzzy_equals
192  ((this->point(13) - this->point(0))*3))
193  return false;
194 
195  return true;
196 }
197 
198 
199 
201 {
202  return THIRD;
203 }
204 
205 
206 
207 unsigned int Tet14::local_side_node(unsigned int side,
208  unsigned int side_node) const
209 {
210  libmesh_assert_less (side, this->n_sides());
211  libmesh_assert_less (side_node, Tet14::nodes_per_side);
212 
213  return Tet14::side_nodes_map[side][side_node];
214 }
215 
216 
217 
218 unsigned int Tet14::local_edge_node(unsigned int edge,
219  unsigned int edge_node) const
220 {
221  libmesh_assert_less (edge, this->n_edges());
222  libmesh_assert_less (edge_node, Tet14::nodes_per_edge);
223 
224  return Tet14::edge_nodes_map[edge][edge_node];
225 }
226 
227 
228 
229 std::unique_ptr<Elem> Tet14::build_side_ptr (const unsigned int i)
230 {
231  return this->simple_build_side_ptr<Tri7, Tet14>(i);
232 }
233 
234 
235 
236 void Tet14::build_side_ptr (std::unique_ptr<Elem> & side,
237  const unsigned int i)
238 {
239  this->simple_build_side_ptr<Tet14>(side, i, TRI7);
240 }
241 
242 
243 
244 std::unique_ptr<Elem> Tet14::build_edge_ptr (const unsigned int i)
245 {
246  return this->simple_build_edge_ptr<Edge3,Tet14>(i);
247 }
248 
249 
250 
251 void Tet14::build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i)
252 {
253  this->simple_build_edge_ptr<Tet14>(edge, i, EDGE3);
254 }
255 
256 
257 
258 void Tet14::connectivity(const unsigned int sc,
259  const IOPackage iop,
260  std::vector<dof_id_type> & conn) const
261 {
263  libmesh_assert_less (sc, 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(8);
271  switch (sc)
272  {
273 
274 
275  // Linear sub-tet 0
276  case 0:
277 
278  conn[0] = this->node_id(0)+1;
279  conn[1] = this->node_id(4)+1;
280  conn[2] = this->node_id(6)+1;
281  conn[3] = this->node_id(6)+1;
282  conn[4] = this->node_id(7)+1;
283  conn[5] = this->node_id(7)+1;
284  conn[6] = this->node_id(7)+1;
285  conn[7] = this->node_id(7)+1;
286 
287  return;
288 
289  // Linear sub-tet 1
290  case 1:
291 
292  conn[0] = this->node_id(4)+1;
293  conn[1] = this->node_id(1)+1;
294  conn[2] = this->node_id(5)+1;
295  conn[3] = this->node_id(5)+1;
296  conn[4] = this->node_id(8)+1;
297  conn[5] = this->node_id(8)+1;
298  conn[6] = this->node_id(8)+1;
299  conn[7] = this->node_id(8)+1;
300 
301  return;
302 
303  // Linear sub-tet 2
304  case 2:
305 
306  conn[0] = this->node_id(5)+1;
307  conn[1] = this->node_id(2)+1;
308  conn[2] = this->node_id(6)+1;
309  conn[3] = this->node_id(6)+1;
310  conn[4] = this->node_id(9)+1;
311  conn[5] = this->node_id(9)+1;
312  conn[6] = this->node_id(9)+1;
313  conn[7] = this->node_id(9)+1;
314 
315  return;
316 
317  // Linear sub-tet 3
318  case 3:
319 
320  conn[0] = this->node_id(7)+1;
321  conn[1] = this->node_id(8)+1;
322  conn[2] = this->node_id(9)+1;
323  conn[3] = this->node_id(9)+1;
324  conn[4] = this->node_id(3)+1;
325  conn[5] = this->node_id(3)+1;
326  conn[6] = this->node_id(3)+1;
327  conn[7] = this->node_id(3)+1;
328 
329  return;
330 
331  // Linear sub-tet 4
332  case 4:
333 
334  conn[0] = this->node_id(4)+1;
335  conn[1] = this->node_id(8)+1;
336  conn[2] = this->node_id(6)+1;
337  conn[3] = this->node_id(6)+1;
338  conn[4] = this->node_id(7)+1;
339  conn[5] = this->node_id(7)+1;
340  conn[6] = this->node_id(7)+1;
341  conn[7] = this->node_id(7)+1;
342 
343  return;
344 
345  // Linear sub-tet 5
346  case 5:
347 
348  conn[0] = this->node_id(4)+1;
349  conn[1] = this->node_id(5)+1;
350  conn[2] = this->node_id(6)+1;
351  conn[3] = this->node_id(6)+1;
352  conn[4] = this->node_id(8)+1;
353  conn[5] = this->node_id(8)+1;
354  conn[6] = this->node_id(8)+1;
355  conn[7] = this->node_id(8)+1;
356 
357  return;
358 
359  // Linear sub-tet 6
360  case 6:
361 
362  conn[0] = this->node_id(5)+1;
363  conn[1] = this->node_id(9)+1;
364  conn[2] = this->node_id(6)+1;
365  conn[3] = this->node_id(6)+1;
366  conn[4] = this->node_id(8)+1;
367  conn[5] = this->node_id(8)+1;
368  conn[6] = this->node_id(8)+1;
369  conn[7] = this->node_id(8)+1;
370 
371  return;
372 
373  // Linear sub-tet 7
374  case 7:
375 
376  conn[0] = this->node_id(7)+1;
377  conn[1] = this->node_id(6)+1;
378  conn[2] = this->node_id(9)+1;
379  conn[3] = this->node_id(9)+1;
380  conn[4] = this->node_id(8)+1;
381  conn[5] = this->node_id(8)+1;
382  conn[6] = this->node_id(8)+1;
383  conn[7] = this->node_id(8)+1;
384 
385  return;
386 
387 
388  default:
389  libmesh_error_msg("Invalid sc = " << sc);
390  }
391  }
392 
393  case VTK:
394  {
395  // VTK has vtkHigherOrderTetra which might have the same
396  // connectivity as our Tet14, but this has not been tested
397  // yet.
398  libmesh_experimental();
399  conn.resize(Tet14::num_nodes);
400  for (auto i : index_range(conn))
401  conn[i] = this->node_id(i);
402  return;
403  }
404 
405  default:
406  libmesh_error_msg("Unsupported IO package " << iop);
407  }
408 }
409 
410 
411 
412 unsigned int Tet14::n_second_order_adjacent_vertices (const unsigned int n) const
413 {
414  switch (n)
415  {
416  case 4:
417  case 5:
418  case 6:
419  case 7:
420  case 8:
421  case 9:
422  return 2;
423 
424  case 10:
425  case 11:
426  case 12:
427  case 13:
428  return 3;
429 
430  default:
431  libmesh_error_msg("Invalid n = " << n);
432  }
433 }
434 
435 
436 
437 
438 const unsigned short int Tet14::_second_order_vertex_child_number[14] =
439  {
440  99,99,99,99, // Vertices
441  0,1,0,0,1,2, // Edges
442  5,4,6,7 // Faces
443  };
444 
445 
446 
447 const unsigned short int Tet14::_second_order_vertex_child_index[14] =
448  {
449  99,99,99,99, // Vertices
450  1,2,2,3,3,3, // Edges
451  10,13,12,13 // Faces
452  };
453 
454 
455 
456 std::pair<unsigned short int, unsigned short int>
457 Tet14::second_order_child_vertex (const unsigned int n) const
458 {
459  libmesh_assert_greater_equal (n, this->n_vertices());
460  libmesh_assert_less (n, this->n_nodes());
461  return std::pair<unsigned short int, unsigned short int>
464 }
465 
466 
467 
468 unsigned short int Tet14::second_order_adjacent_vertex (const unsigned int n,
469  const unsigned int v) const
470 {
471  libmesh_assert_greater_equal (n, this->n_vertices());
472  libmesh_assert_less (n, this->n_nodes());
473  libmesh_assert_less (v, 3);
474  libmesh_assert (n > 9 || v < 2); // Only face nodes have multiple adjacencies
475  return _second_order_adjacent_vertices[n-this->n_vertices()][v];
476 }
477 
478 
479 
480 const unsigned short int Tet14::_second_order_adjacent_vertices[10][3] =
481  {
482  {0, 1, 99}, // vertices adjacent to node 4
483  {1, 2, 99}, // vertices adjacent to node 5
484  {0, 2, 99}, // vertices adjacent to node 6
485  {0, 3, 99}, // vertices adjacent to node 7
486  {1, 3, 99}, // vertices adjacent to node 8
487  {2, 3, 99}, // vertices adjacent to node 9
488  {0, 1, 2}, // vertices adjacent to node 10
489  {0, 1, 3}, // vertices adjacent to node 11
490  {1, 2, 3}, // vertices adjacent to node 12
491  {0, 2, 3}, // vertices adjacent to node 13
492  };
493 
494 
495 
496 
497 
498 #ifdef LIBMESH_ENABLE_AMR
499 
501  {
502  // embedding matrix for child 0
503  {
504  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
505  { 1., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 0
506  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 1
507  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 2
508  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0.}, // 3
509  { 0.375,-0.125, 0., 0., 0.75, 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 4
510  {.09375,-.03125,-.03125, 0., 0.125,-0.125, 0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 5
511  { 0.375, 0.,-0.125, 0., 0., 0., 0.75, 0., 0., 0., 0., 0., 0., 0.}, // 6
512  { 0.375, 0., 0.,-0.125, 0., 0., 0., 0.75, 0., 0., 0., 0., 0., 0.}, // 7
513  {.09375,-.03125, 0.,-.03125,0.125, 0., 0., 0.125,-0.125, 0., 0.,.84375, 0., 0.}, // 8
514  {.09375, 0.,-.03125,-.03125, 0., 0., 0.125, 0.125, 0.,-0.125, 0., 0., 0.,.84375}, // 9
515  { 5/r18,-1/r18,-1/r18, 0., 4/r18,-2/r18, 4/r18, 0., 0., 0., 0.5, 0., 0., 0.}, // 10
516  { 5/r18,-1/r18, 0.,-1/r18, 4/r18, 0., 0., 4/r18,-2/r18, 0., 0., 0.5, 0., 0.}, // 11
517  { 0.125,-1/r72,-1/r72,-1/r72, 0.,-2/r18, 0., 0.,-2/r18,-2/r18, 0.375, 0.375, 0.125, 0.375}, // 12
518  { 5/r18, 0.,-1/r18,-1/r18, 0., 0., 4/r18, 4/r18, 0.,-2/r18, 0., 0., 0., 0.5} // 13
519  },
520 
521  // embedding matrix for child 1?
522  {
523  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
524  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 0
525  { 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 1
526  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0.}, // 2
527  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 3
528  {-0.125, 0.375, 0., 0., 0.75, 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 4
529  { 0., 0.375,-0.125, 0., 0., 0.75, 0., 0., 0., 0., 0., 0., 0., 0.}, // 5
530  {-.03125,.09375,-.03125, 0., 0.125, 0.125,-0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 6
531  {-.03125,.09375, 0.,-.03125,0.125, 0., 0.,-0.125, 0.125, 0., 0.,.84375, 0., 0.}, // 7
532  { 0., 0.375, 0.,-0.125, 0., 0., 0., 0., 0.75, 0., 0., 0., 0., 0.}, // 8
533  { 0.,.09375,-0.03125,-0.03125,0., 0.125, 0., 0., 0.125,-0.125, 0., 0.,.84375, 0.}, // 9
534  {-1/r18, 5/r18,-1/r18, 0., 4/r18, 4/r18,-2/r18, 0., 0., 0., 0.5, 0., 0., 0.}, // 10
535  {-1/r18, 5/r18, 0.,-1/r18, 4/r18, 0., 0.,-2/r18, 4/r18, 0., 0., 0.5, 0., 0.}, // 11
536  { 0., 5/r18,-1/r18,-1/r18, 0., 4/r18, 0., 0., 4/r18,-2/r18, 0., 0., 0.5, 0.}, // 12
537  {-1/r72, 0.125,-1/r72,-1/r72, 0., 0.,-2/r18,-2/r18, 0.,-2/r18, 0.375, 0.375, 0.375, 0.125} // 13
538  },
539 
540  // embedding matrix for child 2
541  {
542  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
543  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 0
544  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0.}, // 1
545  { 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 2
546  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 3
547  {-.03125,-.03125,.09375, 0.,-0.125, 0.125, 0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 4
548  { 0.,-0.125, 0.375, 0., 0., 0.75, 0., 0., 0., 0., 0., 0., 0., 0.}, // 5
549  {-0.125, 0., 0.375, 0., 0., 0., 0.75, 0., 0., 0., 0., 0., 0., 0.}, // 6
550  {-.03125, 0.,.09375,-.03125, 0., 0., 0.125,-0.125, 0., 0.125, 0., 0., 0.,.84375}, // 7
551  { 0.,-.03125,.09375,-.03125, 0., 0.125, 0., 0.,-0.125, 0.125, 0., 0.,.84375, 0.}, // 8
552  { 0., 0., 0.375,-0.125, 0., 0., 0., 0., 0., 0.75, 0., 0., 0., 0.}, // 9
553  {-1/r18,-1/r18, 5/r18, 0.,-2/r18, 4/r18, 4/r18, 0., 0., 0., 0.5, 0., 0., 0.}, // 10
554  {-1/r72,-1/r72, 0.125,-1/r72,-2/r18, 0., 0.,-2/r18,-2/r18, 0., 0.375, 0.125, 0.375, 0.375}, // 11
555  { 0.,-1/r18, 5/r18,-1/r18, 0., 4/r18, 0., 0.,-2/r18, 4/r18, 0., 0., 0.5, 0.}, // 12
556  {-1/r18, 0., 5/r18,-1/r18, 0., 0., 4/r18,-2/r18, 0., 4/r18, 0., 0., 0., 0.5} // 13
557  },
558 
559  // embedding matrix for child 3
560  {
561  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
562  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0.}, // 0
563  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 1
564  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 2
565  { 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 3
566  {-.03125,-.03125, 0.,.09375,-0.125, 0., 0., 0.125, 0.125, 0., 0.,.84375, 0., 0.}, // 4
567  { 0.,-.03125,-.03125,.09375, 0.,-0.125, 0., 0., 0.125, 0.125, 0., 0.,.84375, 0.}, // 5
568  {-.03125, 0.,-.03125,.09375, 0., 0.,-0.125, 0.125, 0., 0.125, 0., 0., 0.,.84375}, // 6
569  {-0.125, 0., 0., 0.375, 0., 0., 0., 0.75, 0., 0., 0., 0., 0., 0.}, // 7
570  { 0.,-0.125, 0., 0.375, 0., 0., 0., 0., 0.75, 0., 0., 0., 0., 0.}, // 8
571  { 0., 0.,-0.125, 0.375, 0., 0., 0., 0., 0., 0.75, 0., 0., 0., 0.}, // 9
572  {-1/r72,-1/r72,-1/r72, 0.125,-2/r18,-2/r18,-2/r18, 0., 0., 0., 0.125, 0.375, 0.375, 0.375}, // 10
573  {-1/r18,-1/r18, 0., 5/r18,-2/r18, 0., 0., 4/r18, 4/r18, 0., 0., 0.5, 0., 0.}, // 11
574  { 0.,-1/r18,-1/r18, 5/r18, 0.,-2/r18, 0., 0., 4/r18, 4/r18, 0., 0., 0.5, 0.}, // 12
575  {-1/r18, 0.,-1/r18, 5/r18, 0., 0.,-2/r18, 4/r18, 0., 4/r18, 0., 0., 0., 0.5} // 13
576  },
577 
578  // embedding matrix for child 4
579  {
580  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
581  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 0
582  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 1
583  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 2
584  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0.}, // 3
585  {-.03125,.09375, 0.,-.03125,0.125, 0., 0.,-0.125, 0.125, 0., 0.,.84375, 0., 0.}, // 4
586  { 1/r64, 1/r64, 1/r64, 1/r64,-0.125,-0.125,-0.125,-0.125,-0.125,-0.125,27/r64,27/r64,27/r64,27/r64}, // 5
587  {.09375,-.03125,-.03125, 0., 0.125,-0.125, 0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 6
588  {.09375,-.03125, 0.,-.03125,0.125, 0., 0., 0.125,-0.125, 0., 0.,.84375, 0., 0.}, // 7
589  {-.03125,-.03125, 0.,.09375,-0.125, 0., 0., 0.125, 0.125, 0., 0.,.84375, 0., 0.}, // 8
590  {.09375, 0.,-.03125,-.03125, 0., 0., 0.125, 0.125, 0.,-0.125, 0., 0., 0.,.84375}, // 9
591  { 2/r72, 2/r72, 0., 0., 0.,-2/r18,-2/r18,-2/r18,-2/r18,-2/r18, 0.5, 0.5, 0.25, 0.25}, // 10
592  { 0., 0, 0., 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 11
593  { 2/r72, 0., 0., 2/r72,-2/r18,-2/r18,-2/r18, 0.,-2/r18,-2/r18, 0.25, 0.5, 0.25, 0.5}, // 12
594  { 0.125,-1/r72,-1/r72,-1/r72, 0.,-2/r18, 0., 0.,-2/r18,-2/r18, 0.375, 0.375, 0.125, 0.375} // 13
595  },
596 
597  // embedding matrix for child 5?
598  {
599  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
600  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 0
601  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0.}, // 1
602  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 2
603  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 3
604  {-.03125,.09375,-.03125, 0., 0.125, 0.125,-0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 4
605  {-.03125,-.03125,.09375, 0.,-0.125, 0.125, 0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 5
606  {.09375,-.03125,-.03125, 0., 0.125,-0.125, 0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 6
607  {-.03125,.09375, 0.,-.03125,0.125, 0., 0.,-0.125, 0.125, 0., 0.,.84375, 0., 0.}, // 7
608  { 0.,.09375,-.03125,-.03125, 0., 0.125, 0., 0., 0.125,-0.125, 0., 0.,.84375, 0.}, // 8
609  { 1/r64, 1/r64, 1/r64, 1/r64,-0.125,-0.125,-0.125,-0.125,-0.125,-0.125,27/r64,27/r64,27/r64,27/r64}, // 9
610  { 0., 0, 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 10
611  {-1/r72, 0.125,-1/r72,-1/r72, 0., 0.,-2/r18,-2/r18, 0.,-2/r18, 0.375, 0.375, 0.375, 0.125}, // 11
612  { 0., 2/r72, 2/r72, 0.,-2/r18, 0.,-2/r18,-2/r18,-2/r18,-2/r18, 0.5, 0.25, 0.5, 0.25}, // 12
613  { 2/r72, 2/r72, 0., 0., 0.,-2/r18,-2/r18,-2/r18,-2/r18,-2/r18, 0.5, 0.5, 0.25, 0.25} // 13
614  },
615 
616  // embedding matrix for child 6?
617  {
618  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
619  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 0
620  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0., 0.}, // 1
621  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 2
622  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 3
623  {-.03125,-.03125,.09375, 0.,-0.125, 0.125, 0.125, 0., 0., 0.,.84375, 0., 0., 0.}, // 4
624  { 0.,-.03125,.09375,-.03125, 0., 0.125, 0., 0.,-0.125, 0.125, 0., 0.,.84375, 0.}, // 5
625  {-.03125, 0.,.09375,-.03125, 0., 0., 0.125,-0.125, 0., 0.125, 0., 0., 0.,.84375}, // 6
626  { 1/r64, 1/r64, 1/r64, 1/r64,-0.125,-0.125,-0.125,-0.125,-0.125,-0.125,27/r64,27/r64,27/r64,27/r64}, // 7
627  { 0.,.09375,-.03125,-.03125, 0., 0.125, 0., 0., 0.125,-0.125, 0., 0.,.84375, 0.}, // 8
628  { 0.,-.03125,-.03125,.09375, 0.,-0.125, 0., 0., 0.125, 0.125, 0., 0.,.84375, 0.}, // 9
629  {-1/r72,-1/r72, 0.125,-1/r72,-2/r18, 0., 0.,-2/r18,-2/r18, 0., 0.375, 0.125, 0.375, 0.375}, // 10
630  { 0., 2/r72, 2/r72, 0.,-2/r18, 0.,-2/r18,-2/r18,-2/r18,-2/r18, 0.5, 0.25, 0.5, 0.25}, // 11
631  { 0., 0, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 12
632  { 0., 0., 2/r72, 2/r72,-2/r18,-2/r18,-2/r18,-2/r18,-2/r18, 0., 0.25, 0.25, 0.5, 0.5} // 13
633  },
634 
635  // embedding matrix for child 7?
636  {
637  // 0 1 2 3 4 5 6 7 8 9 10 11 12 13
638  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 0
639  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 1
640  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 2
641  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0., 0., 0., 0., 0.}, // 3
642  { 1/r64, 1/r64, 1/r64, 1/r64,-0.125,-0.125,-0.125,-0.125,-0.125,-0.125,27/r64,27/r64,27/r64,27/r64}, // 4
643  { 0.,-.03125,-.03125,.09375, 0.,-0.125, 0., 0., 0.125, 0.125, 0., 0.,.84375, 0.}, // 5
644  {-.03125, 0.,.09375,-.03125, 0., 0., 0.125,-0.125, 0., 0.125, 0., 0., 0.,.84375}, // 6
645  {.09375, 0.,-.03125,-.03125, 0., 0., 0.125, 0.125, 0.,-0.125, 0., 0., 0.,.84375}, // 7
646  {-.03125,-.03125, 0.,.09375,-0.125, 0., 0., 0.125, 0.125, 0., 0.,.84375, 0., 0.}, // 8
647  {-.03125, 0.,-.03125,.09375, 0., 0.,-0.125, 0.125, 0., 0.125, 0., 0., 0.,.84375}, // 9
648  { 0., 0., 2/r72, 2/r72,-2/r18,-2/r18,-2/r18,-2/r18,-2/r18, 0., 0.25, 0.25, 0.5, 0.5}, // 10
649  { 2/r72, 0., 0., 2/r72,-2/r18,-2/r18,-2/r18, 0.,-2/r18,-2/r18, 0.25, 0.5, 0.25, 0.5}, // 11
650  {-1/r72,-1/r72,-1/r72, 0.125,-2/r18,-2/r18,-2/r18, 0., 0., 0., 0.125, 0.375, 0.375, 0.375}, // 12
651  { 0., 0, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 13
652  }
653  };
654 
655 
656 
657 Real Tet14::embedding_matrix (const unsigned int i,
658  const unsigned int j,
659  const unsigned int k) const
660 {
661  // Choose an optimal diagonal, if one has not already been selected
662  this->choose_diagonal();
663 
664  // Permuted j and k indices
665  unsigned int
666  jp=j,
667  kp=k;
668 
669  if ((i>3) && (this->_diagonal_selection!=DIAG_02_13))
670  {
671  // Just the enum value cast to an unsigned int...
672  const unsigned ds = static_cast<unsigned int>(this->_diagonal_selection); // == 1 or 2
673 
674  // Instead of doing a lot of arithmetic, use these
675  // straightforward arrays for the permutations. Note that 3 ->
676  // 3 and 10 -> 10, and the first array consists of "forward"
677  // permutations of the sets {0,1,2}, {4,5,6}, {7,8,9}, and
678  // {11,12,13} while the second array consists of "reverse"
679  // permutations of the same sets.
680 
681  const unsigned int perms[2][14] =
682  {
683  {1, 2, 0, 3, 5, 6, 4, 8, 9, 7, 10, 12, 13, 11},
684  {2, 0, 1, 3, 6, 4, 5, 9, 7, 8, 10, 13, 11, 12}
685  };
686 
687  // Permute j
688  jp = perms[ds-1][j];
689  // if (jp<3)
690  // jp = (jp+ds)%3;
691  // else if (jp>3)
692  // jp = (jp-1+ds)%3 + 1 + 3*((jp-1)/3);
693 
694  // Permute k
695  kp = perms[ds-1][k];
696  // if (kp<3)
697  // kp = (kp+ds)%3;
698  // else if (kp>3)
699  // kp = (kp-1+ds)%3 + 1 + 3*((kp-1)/3);
700  }
701 
702  // Debugging:
703  // libMesh::err << "Selected diagonal " << _diagonal_selection << std::endl;
704  // libMesh::err << "j=" << j << std::endl;
705  // libMesh::err << "k=" << k << std::endl;
706  // libMesh::err << "jp=" << jp << std::endl;
707  // libMesh::err << "kp=" << kp << std::endl;
708 
709  // Call embedding matrix with permuted indices
710  return this->_embedding_matrix[i][jp][kp];
711 }
712 
713 
714 const std::vector<std::pair<unsigned char, unsigned char>> &
716  unsigned int n) const
717 {
718  // Choose an optimal diagonal, if one has not already been selected
719  this->choose_diagonal();
720 
721  // Just the enum value cast to an unsigned int...
722  const unsigned ds = static_cast<unsigned int>(this->_diagonal_selection);
723 
724  return Tet14::_parent_bracketing_nodes[ds][c][n];
725 }
726 
727 
728 const std::vector<std::pair<unsigned char, unsigned char>>
730 {
731  // DIAG_02_13
732  {
733  // Node 0, 1, 2, 3, 4, 5, 6, 7, 8, 9,
734  // 10, 11, 12, 13
735  // Child 0
736  { {},{{0,1}},{{0,2}},{{0,3}},{{0,4}},{{4,6}},{{0,6}},{{0,7}},{{4,7}},{{6,7}},{{0,10}},{{0,11}},{{0,12}},{{0,13}} },
737  // Child 1
738  { {{0,1}}, {},{{1,2}},{{1,3}},{{1,4}},{{1,5}},{{4,5}},{{4,8}},{{1,8}},{{5,8}},{{1,10}},{{1,11}},{{1,12}},{{1,13}} },
739  // Child 2
740  { {{0,2}},{{1,2}}, {},{{2,3}},{{5,6}},{{2,5}},{{2,6}},{{6,9}},{{5,9}},{{2,9}},{{2,10}},{{2,11}},{{2,12}},{{2,13}} },
741  // Child 3
742  { {{0,3}},{{1,3}},{{2,3}}, {},{{7,8}},{{8,9}},{{7,9}},{{3,7}},{{3,8}},{{3,9}},{{3,10}},{{3,11}},{{3,12}},{{3,13}} },
743  // Child 4
744  { {{0,1}},{{1,3}},{{0,2}},{{0,3}},{{4,8}},{{6,8}},{{4,6}},{{4,7}},{{7,8}},{{6,7}},
745  {{10,11}},{{0,8},{1,7},{3,4}}, {{11,13}}, {{0,12}} },
746  // Child 5
747  { {{0,1}},{{1,2}},{{0,2}},{{1,3}},{{4,5}},{{5,6}},{{4,6}},{{4,8}},{{5,8}},{{6,8}},
748  {{0,5},{1,6},{2,4}}, {{1,13}}, {{10,12}}, {{10,11}} },
749  // Child 6
750  { {{0,2}},{{1,2}},{{2,3}},{{1,3}},{{5,6}},{{5,9}},{{6,9}},{{6,8}},{{5,8}},{{8,9}},
751  {{2,11}}, {{10,12}},{{1,9},{2,8},{3,5}}, {{12,13}} },
752  // Child 7
753  { {{0,2}},{{1,3}},{{2,3}},{{0,3}},{{6,8}},{{8,9}},{{6,9}},{{6,7}},{{7,8}},{{7,9}},
754  {{12,13}}, {{11,13}}, {{3,10}},{{0,9},{2,7},{3,6}} }
755  },
756  // DIAG_03_12
757  {
758  // Node 0, 1, 2, 3, 4, 5, 6, 7, 8, 9,
759  // 10, 11, 12, 13
760  // Child 0
761  { {},{{0,1}},{{0,2}},{{0,3}},{{0,4}},{{4,6}},{{0,6}},{{0,7}},{{4,7}},{{6,7}},{{0,10}},{{0,11}},{{0,12}},{{0,13}} },
762  // Child 1
763  { {{0,1}}, {},{{1,2}},{{1,3}},{{1,4}},{{1,5}},{{4,5}},{{4,8}},{{1,8}},{{5,8}},{{1,10}},{{1,11}},{{1,12}},{{1,13}} },
764  // Child 2
765  { {{0,2}},{{1,2}}, {},{{2,3}},{{5,6}},{{2,5}},{{2,6}},{{6,9}},{{5,9}},{{2,9}},{{2,10}},{{2,11}},{{2,12}},{{2,13}} },
766  // Child 3
767  { {{0,3}},{{1,3}},{{2,3}}, {},{{7,8}},{{8,9}},{{7,9}},{{3,7}},{{3,8}},{{3,9}},{{3,10}},{{3,11}},{{3,12}},{{3,13}} },
768  // Child 4
769  { {{0,3}},{{1,2}},{{0,2}},{{2,3}},{{5,7}},{{5,6}},{{6,7}},{{7,9}},{{5,9}},{{6,9}},
770  {{10,13}}, {{12,13}}, {{2,11}},{{0,9},{2,7},{3,6}} },
771  // Child 5
772  { {{0,1}},{{1,2}},{{0,2}},{{0,3}},{{4,5}},{{5,6}},{{4,6}},{{4,7}},{{5,7}},{{6,7}},
773  {{0,5},{1,6},{2,4}}, {{10,11}}, {{10,13}}, {{0,12}} },
774  // Child 6
775  { {{0,1}},{{1,3}},{{1,2}},{{0,3}},{{4,8}},{{5,8}},{{4,5}},{{4,7}},{{7,8}},{{5,7}},
776  {{1,13}},{{0,8},{1,7},{3,4}}, {{11,12}}, {{10,11}} },
777  // Child 7
778  { {{0,3}},{{1,3}},{{1,2}},{{2,3}},{{7,8}},{{5,8}},{{5,7}},{{7,9}},{{8,9}},{{5,9}},
779  {{11,12}}, {{3,10}},{{1,9},{2,8},{3,5}}, {{12,13}} }
780  },
781  // DIAG_01_23
782  {
783  // Node 0, 1, 2, 3, 4, 5, 6, 7, 8, 9,
784  // 10, 11, 12, 13
785  // Child 0
786  { {},{{0,1}},{{0,2}},{{0,3}},{{0,4}},{{4,6}},{{0,6}},{{0,7}},{{4,7}},{{6,7}},{{0,10}},{{0,11}},{{0,12}},{{0,13}} },
787  // Child 1
788  { {{0,1}}, {},{{1,2}},{{1,3}},{{1,4}},{{1,5}},{{4,5}},{{4,8}},{{1,8}},{{5,8}},{{1,10}},{{1,11}},{{1,12}},{{1,13}} },
789  // Child 2
790  { {{0,2}},{{1,2}}, {},{{2,3}},{{5,6}},{{2,5}},{{2,6}},{{6,9}},{{5,9}},{{2,9}},{{2,10}},{{2,11}},{{2,12}},{{2,13}} },
791  // Child 3
792  { {{0,3}},{{1,3}},{{2,3}}, {},{{7,8}},{{8,9}},{{7,9}},{{3,7}},{{3,8}},{{3,9}},{{3,10}},{{3,11}},{{3,12}},{{3,13}} },
793  // Child 4
794  { {{0,1}},{{1,2}},{{2,3}},{{1,3}},{{4,5}},{{5,9}},{{4,9}},{{4,8}},{{5,8}},{{8,9}},
795  {{10,12}}, {{1,13}},{{1,9},{2,8},{3,5}}, {{11,12}} },
796  // Child 5
797  { {{0,1}},{{1,2}},{{0,2}},{{2,3}},{{4,5}},{{5,6}},{{4,6}},{{4,9}},{{5,9}},{{6,9}},
798  {{0,5},{1,6},{2,4}}, {{10,12}}, {{2,11}}, {{10,13}} },
799  // Child 6
800  { {{0,3}},{{0,1}},{{0,2}},{{2,3}},{{4,7}},{{4,6}},{{6,7}},{{7,9}},{{4,9}},{{6,9}},
801  {{0,12}}, {{11,13}}, {{10,13}},{{0,9},{2,7},{3,6}} },
802  // Child 7
803  { {{0,3}},{{0,1}},{{2,3}},{{1,3}},{{4,7}},{{4,9}},{{7,9}},{{7,8}},{{4,8}},{{8,9}},
804  {{11,13}},{{0,8},{1,7},{3,4}}, {{11,12}}, {{3,10}} }
805  }
806 };
807 
808 #endif // #ifdef LIBMESH_ENABLE_AMR
809 
810 
811 
812 void Tet14::permute(unsigned int perm_num)
813 {
814  libmesh_assert_less (perm_num, 12);
815 
816  const unsigned int side = perm_num % 4;
817  const unsigned int rotate = perm_num / 4;
818 
819  for (unsigned int i = 0; i != rotate; ++i)
820  {
821  swap3nodes(0,1,2);
822  swap3nodes(4,5,6);
823  swap3nodes(7,8,9);
824  swap3nodes(11,12,13);
825  swap3neighbors(1,2,3);
826  }
827 
828  switch (side) {
829  case 0:
830  break;
831  case 1:
832  swap3nodes(0,2,3);
833  swap3nodes(4,5,8);
834  swap3nodes(6,9,7);
835  swap3nodes(10,12,11);
836  swap3neighbors(0,2,1);
837  break;
838  case 2:
839  swap3nodes(2,0,3);
840  swap3nodes(5,4,8);
841  swap3nodes(6,7,9);
842  swap3nodes(10,11,12);
843  swap3neighbors(0,1,2);
844  break;
845  case 3:
846  swap3nodes(2,1,3);
847  swap3nodes(5,8,9);
848  swap3nodes(6,4,7);
849  swap3nodes(10,11,13);
850  swap3neighbors(0,1,3);
851  break;
852  default:
853  libmesh_error();
854  }
855 }
856 
857 
858 void Tet14::flip(BoundaryInfo * boundary_info)
859 {
860  libmesh_assert(boundary_info);
861 
862  swap2nodes(0,2);
863  swap2nodes(4,5);
864  swap2nodes(7,9);
865  swap2nodes(11,12);
866  swap2neighbors(1,2);
867  swap2boundarysides(1,2,boundary_info);
868  swap2boundaryedges(0,1,boundary_info);
869  swap2boundaryedges(3,5,boundary_info);
870 }
871 
872 
873 ElemType Tet14::side_type (const unsigned int libmesh_dbg_var(s)) const
874 {
875  libmesh_assert_less (s, 4);
876  return TRI7;
877 }
878 
879 
880 } // namespace libMesh
ElemType side_type(const unsigned int s) const override final
Definition: cell_tet14.C:873
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
static const Real _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
Definition: cell_tet14.h:272
Node ** _nodes
Pointers to the nodes we are connected to.
Definition: elem.h:2245
virtual Order default_order() const override
Definition: cell_tet14.C:200
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: cell_tet14.C:97
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const override
Definition: cell_tet14.C:468
virtual unsigned int n_nodes() const override
Definition: cell_tet14.h:96
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
This maps the node of the edge to element node numbers.
Definition: cell_tet14.h:235
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
Definition: cell_tet14.C:457
virtual bool is_face(const unsigned int i) const override
Definition: cell_tet14.C:80
virtual unsigned int n_sides() const override final
Definition: cell_tet.h:81
static const int num_sides
Geometric constants for all Tets.
Definition: cell_tet.h:74
static const std::vector< std::pair< unsigned char, unsigned char > > _parent_bracketing_nodes[3][num_children][num_nodes]
Pairs of nodes that bracket child nodes when doing mesh refinement, for each of the three possible di...
Definition: cell_tet14.h:280
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
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 int num_edges
Definition: cell_tet.h:75
static const int nodes_per_edge
Definition: cell_tet14.h:223
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.
Definition: cell_tet14.C:657
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Builds a TRI7 built coincident with face i.
Definition: cell_tet14.C:229
Diagonal _diagonal_selection
The currently-selected diagonal used during refinement.
Definition: cell_tet.h:249
static const unsigned short int _second_order_vertex_child_index[14]
Vector that names the child vertex index for each second order node.
Definition: cell_tet14.h:302
The libMesh namespace provides an interface to certain functionality in the library.
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
Definition: cell_tet14.C:104
void swap3neighbors(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three neighbor_ptrs, "rotating" them.
Definition: elem.h:2133
static const unsigned short int _second_order_adjacent_vertices[10][3]
Matrix that tells which vertices define the location of mid-side or mid-face nodes, indexed by node_num-4.
Definition: cell_tet14.h:292
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: cell_tet14.C:258
virtual unsigned int n_vertices() const override final
Definition: cell_tet.h:86
virtual void permute(unsigned int perm_num) override final
Permutes the element (by swapping node and neighbor pointers) according to the specified index...
Definition: cell_tet14.C:812
void swap3nodes(unsigned int n1, unsigned int n2, unsigned int n3)
Swaps three node_ptrs, "rotating" them.
Definition: elem.h:2124
virtual unsigned int n_sub_elem() const override
Definition: cell_tet14.h:101
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
virtual bool has_affine_map() const override
Definition: cell_tet14.C:153
virtual bool is_edge(const unsigned int i) const override
Definition: cell_tet14.C:73
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: cell_tet14.C:110
virtual unsigned int n_edges() const override final
Definition: cell_tet.h:91
RealTensorValue rotate(MeshBase &mesh, const Real phi, const Real theta=0., const Real psi=0.)
Rotates the mesh in the xy plane.
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:2102
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Definition: cell_tet14.C:207
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE3 built coincident with edge i.
Definition: cell_tet14.C:244
static const int num_nodes
Geometric constants for Tet14.
Definition: cell_tet14.h:221
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
Definition: cell_tet14.C:123
static const int nodes_per_side
Definition: cell_tet14.h:222
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool is_vertex(const unsigned int i) const override
Definition: cell_tet14.C:66
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
Definition: cell_tet14.h:229
static const unsigned short int _second_order_vertex_child_number[14]
Vector that names a child sharing each second order node.
Definition: cell_tet14.h:297
void choose_diagonal() const
Derived classes use this function to select an initial diagonal during refinement.
Definition: cell_tet.C:204
static const int num_children
Definition: cell_tet.h:76
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
Definition: cell_tet14.C:858
virtual unsigned int n_second_order_adjacent_vertices(const unsigned int) const override
Definition: cell_tet14.C:412
virtual const std::vector< std::pair< unsigned char, unsigned char > > & parent_bracketing_nodes(unsigned int c, unsigned int n) const override
Definition: cell_tet14.C:715
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
Definition: cell_tet14.C:218
bool is_child_on_side_helper(const unsigned int c, const unsigned int s, const unsigned int checked_nodes[][3]) const
Called by descendant classes with appropriate data to determine if child c is on side s...
Definition: cell_tet.C:152
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
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 bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: cell_tet14.C:87