libMesh
cell_tet10.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_tet10.h"
21 #include "libmesh/edge_edge3.h"
22 #include "libmesh/face_tri6.h"
23 #include "libmesh/enum_io_package.h"
24 #include "libmesh/enum_order.h"
25 
26 namespace libMesh
27 {
28 
29 
30 
31 // ------------------------------------------------------------
32 // Tet10 class static member initializations
33 const int Tet10::num_nodes;
34 const int Tet10::nodes_per_side;
35 const int Tet10::nodes_per_edge;
36 
38  {
39  {0, 2, 1, 6, 5, 4}, // Side 0
40  {0, 1, 3, 4, 8, 7}, // Side 1
41  {1, 2, 3, 5, 9, 8}, // Side 2
42  {2, 0, 3, 6, 7, 9} // Side 3
43  };
44 
46  {
47  {0, 1, 4}, // Edge 0
48  {1, 2, 5}, // Edge 1
49  {0, 2, 6}, // Edge 2
50  {0, 3, 7}, // Edge 3
51  {1, 3, 8}, // Edge 4
52  {2, 3, 9} // Edge 5
53  };
54 
55 // ------------------------------------------------------------
56 // Tet10 class member functions
57 
58 bool Tet10::is_vertex(const unsigned int i) const
59 {
60  if (i < 4)
61  return true;
62  return false;
63 }
64 
65 bool Tet10::is_edge(const unsigned int i) const
66 {
67  if (i < 4)
68  return false;
69  return true;
70 }
71 
72 bool Tet10::is_face(const unsigned int) const
73 {
74  return false;
75 }
76 
77 bool Tet10::is_node_on_side(const unsigned int n,
78  const unsigned int s) const
79 {
80  libmesh_assert_less (s, n_sides());
81  return std::find(std::begin(side_nodes_map[s]),
82  std::end(side_nodes_map[s]),
83  n) != std::end(side_nodes_map[s]);
84 }
85 
86 std::vector<unsigned>
87 Tet10::nodes_on_side(const unsigned int s) const
88 {
89  libmesh_assert_less(s, n_sides());
90  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
91 }
92 
93 std::vector<unsigned>
94 Tet10::nodes_on_edge(const unsigned int e) const
95 {
96  libmesh_assert_less(e, n_edges());
97  return {std::begin(edge_nodes_map[e]), std::end(edge_nodes_map[e])};
98 }
99 
100 bool Tet10::is_node_on_edge(const unsigned int n,
101  const unsigned int e) const
102 {
103  libmesh_assert_less (e, n_edges());
104  return std::find(std::begin(edge_nodes_map[e]),
105  std::end(edge_nodes_map[e]),
106  n) != std::end(edge_nodes_map[e]);
107 }
108 
109 
110 #ifdef LIBMESH_ENABLE_AMR
111 
112 // This function only works if LIBMESH_ENABLE_AMR...
113 bool Tet10::is_child_on_side(const unsigned int c,
114  const unsigned int s) const
115 {
116  // Table of local IDs for the midege nodes on the side opposite a given node.
117  // See the ASCII art in the header file for this class to confirm this.
118  const unsigned int midedge_nodes_opposite[4][3] =
119  {
120  {5,8,9}, // midedge nodes opposite node 0
121  {6,7,9}, // midedge nodes opposite node 1
122  {4,7,8}, // midedge nodes opposite node 2
123  {4,5,6} // midedge nodes opposite node 3
124  };
125 
126  // Call the base class helper function
127  return Tet::is_child_on_side_helper(c, s, midedge_nodes_opposite);
128 }
129 
130 #else
131 
132 bool Tet10::is_child_on_side(const unsigned int /*c*/,
133  const unsigned int /*s*/) const
134 {
135  libmesh_not_implemented();
136  return false;
137 }
138 
139 #endif //LIBMESH_ENABLE_AMR
140 
141 
142 
144 {
145  // Make sure edges are straight
146  Point v = this->point(1) - this->point(0);
147  if (!v.relative_fuzzy_equals
148  ((this->point(4) - this->point(0))*2, affine_tol))
149  return false;
150  v = this->point(2) - this->point(1);
151  if (!v.relative_fuzzy_equals
152  ((this->point(5) - this->point(1))*2, affine_tol))
153  return false;
154  v = this->point(2) - this->point(0);
155  if (!v.relative_fuzzy_equals
156  ((this->point(6) - this->point(0))*2, affine_tol))
157  return false;
158  v = this->point(3) - this->point(0);
159  if (!v.relative_fuzzy_equals
160  ((this->point(7) - this->point(0))*2, affine_tol))
161  return false;
162  v = this->point(3) - this->point(1);
163  if (!v.relative_fuzzy_equals
164  ((this->point(8) - this->point(1))*2, affine_tol))
165  return false;
166  v = this->point(3) - this->point(2);
167  if (!v.relative_fuzzy_equals
168  ((this->point(9) - this->point(2))*2, affine_tol))
169  return false;
170  return true;
171 }
172 
173 
174 
176 {
177  return SECOND;
178 }
179 
180 
181 
182 unsigned int Tet10::local_side_node(unsigned int side,
183  unsigned int side_node) const
184 {
185  libmesh_assert_less (side, this->n_sides());
186  libmesh_assert_less (side_node, Tet10::nodes_per_side);
187 
188  return Tet10::side_nodes_map[side][side_node];
189 }
190 
191 
192 
193 unsigned int Tet10::local_edge_node(unsigned int edge,
194  unsigned int edge_node) const
195 {
196  libmesh_assert_less (edge, this->n_edges());
197  libmesh_assert_less (edge_node, Tet10::nodes_per_edge);
198 
199  return Tet10::edge_nodes_map[edge][edge_node];
200 }
201 
202 
203 
204 std::unique_ptr<Elem> Tet10::build_side_ptr (const unsigned int i)
205 {
206  return this->simple_build_side_ptr<Tri6, Tet10>(i);
207 }
208 
209 
210 
211 void Tet10::build_side_ptr (std::unique_ptr<Elem> & side,
212  const unsigned int i)
213 {
214  this->simple_build_side_ptr<Tet10>(side, i, TRI6);
215 }
216 
217 
218 
219 std::unique_ptr<Elem> Tet10::build_edge_ptr (const unsigned int i)
220 {
221  return this->simple_build_edge_ptr<Edge3,Tet10>(i);
222 }
223 
224 
225 
226 void Tet10::build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i)
227 {
228  this->simple_build_edge_ptr<Tet10>(edge, i, EDGE3);
229 }
230 
231 
232 
233 void Tet10::connectivity(const unsigned int sc,
234  const IOPackage iop,
235  std::vector<dof_id_type> & conn) const
236 {
238  libmesh_assert_less (sc, this->n_sub_elem());
239  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
240 
241  switch (iop)
242  {
243  case TECPLOT:
244  {
245  conn.resize(8);
246  switch (sc)
247  {
248 
249 
250  // Linear sub-tet 0
251  case 0:
252 
253  conn[0] = this->node_id(0)+1;
254  conn[1] = this->node_id(4)+1;
255  conn[2] = this->node_id(6)+1;
256  conn[3] = this->node_id(6)+1;
257  conn[4] = this->node_id(7)+1;
258  conn[5] = this->node_id(7)+1;
259  conn[6] = this->node_id(7)+1;
260  conn[7] = this->node_id(7)+1;
261 
262  return;
263 
264  // Linear sub-tet 1
265  case 1:
266 
267  conn[0] = this->node_id(4)+1;
268  conn[1] = this->node_id(1)+1;
269  conn[2] = this->node_id(5)+1;
270  conn[3] = this->node_id(5)+1;
271  conn[4] = this->node_id(8)+1;
272  conn[5] = this->node_id(8)+1;
273  conn[6] = this->node_id(8)+1;
274  conn[7] = this->node_id(8)+1;
275 
276  return;
277 
278  // Linear sub-tet 2
279  case 2:
280 
281  conn[0] = this->node_id(5)+1;
282  conn[1] = this->node_id(2)+1;
283  conn[2] = this->node_id(6)+1;
284  conn[3] = this->node_id(6)+1;
285  conn[4] = this->node_id(9)+1;
286  conn[5] = this->node_id(9)+1;
287  conn[6] = this->node_id(9)+1;
288  conn[7] = this->node_id(9)+1;
289 
290  return;
291 
292  // Linear sub-tet 3
293  case 3:
294 
295  conn[0] = this->node_id(7)+1;
296  conn[1] = this->node_id(8)+1;
297  conn[2] = this->node_id(9)+1;
298  conn[3] = this->node_id(9)+1;
299  conn[4] = this->node_id(3)+1;
300  conn[5] = this->node_id(3)+1;
301  conn[6] = this->node_id(3)+1;
302  conn[7] = this->node_id(3)+1;
303 
304  return;
305 
306  // Linear sub-tet 4
307  case 4:
308 
309  conn[0] = this->node_id(4)+1;
310  conn[1] = this->node_id(8)+1;
311  conn[2] = this->node_id(6)+1;
312  conn[3] = this->node_id(6)+1;
313  conn[4] = this->node_id(7)+1;
314  conn[5] = this->node_id(7)+1;
315  conn[6] = this->node_id(7)+1;
316  conn[7] = this->node_id(7)+1;
317 
318  return;
319 
320  // Linear sub-tet 5
321  case 5:
322 
323  conn[0] = this->node_id(4)+1;
324  conn[1] = this->node_id(5)+1;
325  conn[2] = this->node_id(6)+1;
326  conn[3] = this->node_id(6)+1;
327  conn[4] = this->node_id(8)+1;
328  conn[5] = this->node_id(8)+1;
329  conn[6] = this->node_id(8)+1;
330  conn[7] = this->node_id(8)+1;
331 
332  return;
333 
334  // Linear sub-tet 6
335  case 6:
336 
337  conn[0] = this->node_id(5)+1;
338  conn[1] = this->node_id(9)+1;
339  conn[2] = this->node_id(6)+1;
340  conn[3] = this->node_id(6)+1;
341  conn[4] = this->node_id(8)+1;
342  conn[5] = this->node_id(8)+1;
343  conn[6] = this->node_id(8)+1;
344  conn[7] = this->node_id(8)+1;
345 
346  return;
347 
348  // Linear sub-tet 7
349  case 7:
350 
351  conn[0] = this->node_id(7)+1;
352  conn[1] = this->node_id(6)+1;
353  conn[2] = this->node_id(9)+1;
354  conn[3] = this->node_id(9)+1;
355  conn[4] = this->node_id(8)+1;
356  conn[5] = this->node_id(8)+1;
357  conn[6] = this->node_id(8)+1;
358  conn[7] = this->node_id(8)+1;
359 
360  return;
361 
362 
363  default:
364  libmesh_error_msg("Invalid sc = " << sc);
365  }
366  }
367 
368  case VTK:
369  {
370  // VTK connectivity for VTK_QUADRATIC_TETRA matches libMesh's own.
371  conn.resize(Tet10::num_nodes);
372  for (auto i : index_range(conn))
373  conn[i] = this->node_id(i);
374  return;
375  }
376 
377  default:
378  libmesh_error_msg("Unsupported IO package " << iop);
379  }
380 }
381 
382 
383 
384 const unsigned short int Tet10::_second_order_vertex_child_number[10] =
385  {
386  99,99,99,99, // Vertices
387  0,1,0,0,1,2 // Edges
388  };
389 
390 
391 
392 const unsigned short int Tet10::_second_order_vertex_child_index[10] =
393  {
394  99,99,99,99, // Vertices
395  1,2,2,3,3,3 // Edges
396  };
397 
398 
399 
400 std::pair<unsigned short int, unsigned short int>
401 Tet10::second_order_child_vertex (const unsigned int n) const
402 {
403  libmesh_assert_greater_equal (n, this->n_vertices());
404  libmesh_assert_less (n, this->n_nodes());
405  return std::pair<unsigned short int, unsigned short int>
408 }
409 
410 
411 
412 unsigned short int Tet10::second_order_adjacent_vertex (const unsigned int n,
413  const unsigned int v) const
414 {
415  libmesh_assert_greater_equal (n, this->n_vertices());
416  libmesh_assert_less (n, this->n_nodes());
417  libmesh_assert_less (v, 2);
418  return _second_order_adjacent_vertices[n-this->n_vertices()][v];
419 }
420 
421 
422 
423 const unsigned short int Tet10::_second_order_adjacent_vertices[6][2] =
424  {
425  {0, 1}, // vertices adjacent to node 4
426  {1, 2}, // vertices adjacent to node 5
427  {0, 2}, // vertices adjacent to node 6
428  {0, 3}, // vertices adjacent to node 7
429  {1, 3}, // vertices adjacent to node 8
430  {2, 3} // vertices adjacent to node 9
431  };
432 
433 
434 
435 
436 
437 #ifdef LIBMESH_ENABLE_AMR
438 
440  {
441  // embedding matrix for child 0
442  {
443  // 0 1 2 3 4 5 6 7 8 9
444  { 1., 0., 0., 0., 0., 0., 0., 0., 0., 0.}, // 0
445  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 1
446  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 2
447  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 3
448  { 0.375,-0.125, 0., 0., 0.75, 0., 0., 0., 0., 0.}, // 4
449  { 0.,-0.125,-0.125, 0., 0.5, 0.25, 0.5, 0., 0., 0.}, // 5
450  { 0.375, 0.,-0.125, 0., 0., 0., 0.75, 0., 0., 0.}, // 6
451  { 0.375, 0., 0.,-0.125, 0., 0., 0., 0.75, 0., 0.}, // 7
452  { 0.,-0.125, 0.,-0.125, 0.5, 0., 0., 0.5, 0.25, 0.}, // 8
453  { 0., 0.,-0.125,-0.125, 0., 0., 0.5, 0.5, 0., 0.25} // 9
454  },
455 
456  // embedding matrix for child 1
457  {
458  // 0 1 2 3 4 5 6 7 8 9
459  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 0
460  { 0., 1., 0., 0., 0., 0., 0., 0., 0., 0.}, // 1
461  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 2
462  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 3
463  {-0.125, 0.375, 0., 0., 0.75, 0., 0., 0., 0., 0.}, // 4
464  { 0., 0.375,-0.125, 0., 0., 0.75, 0., 0., 0., 0.}, // 5
465  {-0.125, 0.,-0.125, 0., 0.5, 0.5, 0.25, 0., 0., 0.}, // 6
466  {-0.125, 0., 0.,-0.125, 0.5, 0., 0., 0.25, 0.5, 0.}, // 7
467  { 0., 0.375, 0.,-0.125, 0., 0., 0., 0., 0.75, 0.}, // 8
468  { 0., 0.,-0.125,-0.125, 0., 0.5, 0., 0., 0.5, 0.25} // 9
469  },
470 
471  // embedding matrix for child 2
472  {
473  // 0 1 2 3 4 5 6 7 8 9
474  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 0
475  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 1
476  { 0., 0., 1., 0., 0., 0., 0., 0., 0., 0.}, // 2
477  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 3
478  {-0.125,-0.125, 0., 0., 0.25, 0.5, 0.5, 0., 0., 0.}, // 4
479  { 0.,-0.125, 0.375, 0., 0., 0.75, 0., 0., 0., 0.}, // 5
480  {-0.125, 0., 0.375, 0., 0., 0., 0.75, 0., 0., 0.}, // 6
481  {-0.125, 0., 0.,-0.125, 0., 0., 0.5, 0.25, 0., 0.5}, // 7
482  { 0.,-0.125, 0.,-0.125, 0., 0.5, 0., 0., 0.25, 0.5}, // 8
483  { 0., 0., 0.375,-0.125, 0., 0., 0., 0., 0., 0.75} // 9
484  },
485 
486  // embedding matrix for child 3
487  {
488  // 0 1 2 3 4 5 6 7 8 9
489  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 0
490  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 1
491  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 2
492  { 0., 0., 0., 1., 0., 0., 0., 0., 0., 0.}, // 3
493  {-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5, 0.}, // 4
494  { 0.,-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5}, // 5
495  {-0.125, 0.,-0.125, 0., 0., 0., 0.25, 0.5, 0., 0.5}, // 6
496  {-0.125, 0., 0., 0.375, 0., 0., 0., 0.75, 0., 0.}, // 7
497  { 0.,-0.125, 0., 0.375, 0., 0., 0., 0., 0.75, 0.}, // 8
498  { 0., 0.,-0.125, 0.375, 0., 0., 0., 0., 0., 0.75} // 9
499  },
500 
501  // embedding matrix for child 4
502  {
503  // 0 1 2 3 4 5 6 7 8 9
504  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 0
505  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 1
506  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 2
507  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 3
508  {-0.125, 0., 0.,-0.125, 0.5, 0., 0., 0.25, 0.5, 0.}, // 4
509  {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25}, // 5
510  { 0.,-0.125,-0.125, 0., 0.5, 0.25, 0.5, 0., 0., 0.}, // 6
511  { 0.,-0.125, 0.,-0.125, 0.5, 0., 0., 0.5, 0.25, 0.}, // 7
512  {-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5, 0.}, // 8
513  { 0., 0.,-0.125,-0.125, 0., 0., 0.5, 0.5, 0., 0.25} // 9
514  },
515 
516  // embedding matrix for child 5
517  {
518  // 0 1 2 3 4 5 6 7 8 9
519  { 0., 0., 0., 0., 1., 0., 0., 0., 0., 0.}, // 0
520  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 1
521  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 2
522  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 3
523  {-0.125, 0.,-0.125, 0., 0.5, 0.5, 0.25, 0., 0., 0.}, // 4
524  {-0.125,-0.125, 0., 0., 0.25, 0.5, 0.5, 0., 0., 0.}, // 5
525  { 0.,-0.125,-0.125, 0., 0.5, 0.25, 0.5, 0., 0., 0.}, // 6
526  {-0.125, 0., 0.,-0.125, 0.5, 0., 0., 0.25, 0.5, 0.}, // 7
527  { 0., 0.,-0.125,-0.125, 0., 0.5, 0., 0., 0.5, 0.25}, // 8
528  {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25} // 9
529  },
530 
531  // embedding matrix for child 6
532  {
533  // 0 1 2 3 4 5 6 7 8 9
534  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 0
535  { 0., 0., 0., 0., 0., 1., 0., 0., 0., 0.}, // 1
536  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 2
537  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 3
538  {-0.125,-0.125, 0., 0., 0.25, 0.5, 0.5, 0., 0., 0.}, // 4
539  { 0.,-0.125, 0.,-0.125, 0., 0.5, 0., 0., 0.25, 0.5}, // 5
540  {-0.125, 0., 0.,-0.125, 0., 0., 0.5, 0.25, 0., 0.5}, // 6
541  {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25}, // 7
542  { 0., 0.,-0.125,-0.125, 0., 0.5, 0., 0., 0.5, 0.25}, // 8
543  { 0.,-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5} // 9
544  },
545 
546  // embedding matrix for child 7
547  {
548  // 0 1 2 3 4 5 6 7 8 9
549  { 0., 0., 0., 0., 0., 0., 1., 0., 0., 0.}, // 0
550  { 0., 0., 0., 0., 0., 0., 0., 0., 1., 0.}, // 1
551  { 0., 0., 0., 0., 0., 0., 0., 0., 0., 1.}, // 2
552  { 0., 0., 0., 0., 0., 0., 0., 1., 0., 0.}, // 3
553  {-0.125,-0.125,-0.125,-0.125, 0.25, 0.25, 0.25, 0.25, 0.25, 0.25}, // 4
554  { 0.,-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5}, // 5
555  {-0.125, 0., 0.,-0.125, 0., 0., 0.5, 0.25, 0., 0.5}, // 6
556  { 0., 0.,-0.125,-0.125, 0., 0., 0.5, 0.5, 0., 0.25}, // 7
557  {-0.125,-0.125, 0., 0., 0.25, 0., 0., 0.5, 0.5, 0.}, // 8
558  {-0.125, 0.,-0.125, 0., 0., 0., 0.25, 0.5, 0., 0.5} // 9
559  }
560  };
561 
562 
563 
564 Real Tet10::embedding_matrix (const unsigned int i,
565  const unsigned int j,
566  const unsigned int k) const
567 {
568  // Choose an optimal diagonal, if one has not already been selected
569  this->choose_diagonal();
570 
571  // Permuted j and k indices
572  unsigned int
573  jp=j,
574  kp=k;
575 
576  if ((i>3) && (this->_diagonal_selection!=DIAG_02_13))
577  {
578  // Just the enum value cast to an unsigned int...
579  const unsigned ds = static_cast<unsigned int>(this->_diagonal_selection); // == 1 or 2
580 
581  // Instead of doing a lot of arithmetic, use these
582  // straightforward arrays for the permutations. Note that 3 ->
583  // 3, and the first array consists of "forward" permutations of
584  // the sets {0,1,2}, {4,5,6}, and {7,8,9} while the second array
585  // consists of "reverse" permutations of the same sets.
586  const unsigned int perms[2][10] =
587  {
588  {1, 2, 0, 3, 5, 6, 4, 8, 9, 7},
589  {2, 0, 1, 3, 6, 4, 5, 9, 7, 8}
590  };
591 
592  // Permute j
593  jp = perms[ds-1][j];
594  // if (jp<3)
595  // jp = (jp+ds)%3;
596  // else if (jp>3)
597  // jp = (jp-1+ds)%3 + 1 + 3*((jp-1)/3);
598 
599  // Permute k
600  kp = perms[ds-1][k];
601  // if (kp<3)
602  // kp = (kp+ds)%3;
603  // else if (kp>3)
604  // kp = (kp-1+ds)%3 + 1 + 3*((kp-1)/3);
605  }
606 
607  // Debugging:
608  // libMesh::err << "Selected diagonal " << _diagonal_selection << std::endl;
609  // libMesh::err << "j=" << j << std::endl;
610  // libMesh::err << "k=" << k << std::endl;
611  // libMesh::err << "jp=" << jp << std::endl;
612  // libMesh::err << "kp=" << kp << std::endl;
613 
614  // Call embedding matrix with permuted indices
615  return this->_embedding_matrix[i][jp][kp];
616 }
617 
618 #endif // #ifdef LIBMESH_ENABLE_AMR
619 
620 
621 
623 {
624  // This specialization is good for Lagrange mappings only
625  if (this->mapping_type() != LAGRANGE_MAP)
626  return this->Elem::volume();
627 
628  // Make copies of our points. It makes the subsequent calculations a bit
629  // shorter and avoids dereferencing the same pointer multiple times.
630  Point
631  x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4),
632  x5 = point(5), x6 = point(6), x7 = point(7), x8 = point(8), x9 = point(9);
633 
634  // The constant components of the dx/dxi vector, linear in xi, eta, zeta.
635  // These were copied directly from the output of a Python script.
636  Point dx_dxi[4] =
637  {
638  -3*x0 - x1 + 4*x4, // constant
639  4*x0 - 4*x4 - 4*x7 + 4*x8, // zeta
640  4*x0 - 4*x4 + 4*x5 - 4*x6, // eta
641  4*x0 + 4*x1 - 8*x4 // xi
642  };
643 
644  // The constant components of the dx/deta vector, linear in xi, eta, zeta.
645  // These were copied directly from the output of a Python script.
646  Point dx_deta[4] =
647  {
648  -3*x0 - x2 + 4*x6, // constant
649  4*x0 - 4*x6 - 4*x7 + 4*x9, // zeta
650  4*x0 + 4*x2 - 8*x6, // eta
651  4*x0 - 4*x4 + 4*x5 - 4*x6 // xi
652  };
653 
654  // The constant components of the dx/dzeta vector, linear in xi, eta, zeta.
655  // These were copied directly from the output of a Python script.
656  Point dx_dzeta[4] =
657  {
658  -3*x0 - x3 + 4*x7, // constant
659  4*x0 + 4*x3 - 8*x7, // zeta
660  4*x0 - 4*x6 - 4*x7 + 4*x9, // eta
661  4*x0 - 4*x4 - 4*x7 + 4*x8 // xi
662  };
663 
664  // 2x2x2 conical quadrature rule. Note: there is also a five point
665  // rule for tets with a negative weight which would be cheaper, but
666  // we'll use this one to preclude any possible issues with
667  // cancellation error.
668  const int N = 8;
669  static const Real w[N] =
670  {
671  3.6979856358852914509238091810505e-02_R,
672  1.6027040598476613723156741868689e-02_R,
673  2.1157006454524061178256145400082e-02_R,
674  9.1694299214797439226823542540576e-03_R,
675  3.6979856358852914509238091810505e-02_R,
676  1.6027040598476613723156741868689e-02_R,
677  2.1157006454524061178256145400082e-02_R,
678  9.1694299214797439226823542540576e-03_R
679  };
680 
681  static const Real xi[N] =
682  {
683  1.2251482265544137786674043037115e-01_R,
684  5.4415184401122528879992623629551e-01_R,
685  1.2251482265544137786674043037115e-01_R,
686  5.4415184401122528879992623629551e-01_R,
687  1.2251482265544137786674043037115e-01_R,
688  5.4415184401122528879992623629551e-01_R,
689  1.2251482265544137786674043037115e-01_R,
690  5.4415184401122528879992623629551e-01_R
691  };
692 
693  static const Real eta[N] =
694  {
695  1.3605497680284601717109468420738e-01_R,
696  7.0679724159396903069267439165167e-02_R,
697  5.6593316507280088053551297149570e-01_R,
698  2.9399880063162286589079157179842e-01_R,
699  1.3605497680284601717109468420738e-01_R,
700  7.0679724159396903069267439165167e-02_R,
701  5.6593316507280088053551297149570e-01_R,
702  2.9399880063162286589079157179842e-01_R
703  };
704 
705  static const Real zeta[N] =
706  {
707  1.5668263733681830907933725249176e-01_R,
708  8.1395667014670255076709592007207e-02_R,
709  6.5838687060044409936029672711329e-02_R,
710  3.4202793236766414300604458388142e-02_R,
711  5.8474756320489429588282763292971e-01_R,
712  3.0377276481470755305409673253211e-01_R,
713  2.4571332521171333166171692542182e-01_R,
714  1.2764656212038543100867773351792e-01_R
715  };
716 
717  Real vol = 0.;
718  for (int q=0; q<N; ++q)
719  {
720  // Compute dx_dxi, dx_deta, dx_dzeta at the current quadrature point.
721  Point
722  dx_dxi_q = dx_dxi[0] + zeta[q]*dx_dxi[1] + eta[q]*dx_dxi[2] + xi[q]*dx_dxi[3],
723  dx_deta_q = dx_deta[0] + zeta[q]*dx_deta[1] + eta[q]*dx_deta[2] + xi[q]*dx_deta[3],
724  dx_dzeta_q = dx_dzeta[0] + zeta[q]*dx_dzeta[1] + eta[q]*dx_dzeta[2] + xi[q]*dx_dzeta[3];
725 
726  // Compute scalar triple product, multiply by weight, and accumulate volume.
727  vol += w[q] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q);
728  }
729 
730  return vol;
731 }
732 
733 
734 void Tet10::permute(unsigned int perm_num)
735 {
736  libmesh_assert_less (perm_num, 12);
737 
738  const unsigned int side = perm_num % 4;
739  const unsigned int rotate = perm_num / 4;
740 
741  for (unsigned int i = 0; i != rotate; ++i)
742  {
743  swap3nodes(0,1,2);
744  swap3nodes(4,5,6);
745  swap3nodes(7,8,9);
746  swap3neighbors(1,2,3);
747  }
748 
749  switch (side) {
750  case 0:
751  break;
752  case 1:
753  swap3nodes(0,2,3);
754  swap3nodes(4,5,8);
755  swap3nodes(6,9,7);
756  swap3neighbors(0,2,1);
757  break;
758  case 2:
759  swap3nodes(2,0,3);
760  swap3nodes(5,4,8);
761  swap3nodes(6,7,9);
762  swap3neighbors(0,1,2);
763  break;
764  case 3:
765  swap3nodes(2,1,3);
766  swap3nodes(5,8,9);
767  swap3nodes(6,4,7);
768  swap3neighbors(0,1,3);
769  break;
770  default:
771  libmesh_error();
772  }
773 }
774 
775 
776 void Tet10::flip(BoundaryInfo * boundary_info)
777 {
778  libmesh_assert(boundary_info);
779 
780  swap2nodes(0,2);
781  swap2nodes(4,5);
782  swap2nodes(7,9);
783  swap2neighbors(1,2);
784  swap2boundarysides(1,2,boundary_info);
785  swap2boundaryedges(0,1,boundary_info);
786  swap2boundaryedges(3,5,boundary_info);
787 }
788 
789 
790 ElemType Tet10::side_type (const unsigned int libmesh_dbg_var(s)) const
791 {
792  libmesh_assert_less (s, 4);
793  return TRI6;
794 }
795 
796 
797 } // 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
virtual Order default_order() const override
Definition: cell_tet10.C:175
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
Node ** _nodes
Pointers to the nodes we are connected to.
Definition: elem.h:2245
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE3 built coincident with edge i.
Definition: cell_tet10.C:219
virtual unsigned int n_nodes() const override
Definition: cell_tet10.h:90
virtual bool has_affine_map() const override
Definition: cell_tet10.C:143
virtual unsigned int n_sides() const override final
Definition: cell_tet.h:81
static const unsigned short int _second_order_adjacent_vertices[6][2]
Matrix that tells which vertices define the location of mid-side (or second-order) nodes...
Definition: cell_tet10.h:277
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_tet10.h:230
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const override
Definition: cell_tet10.C:412
static const int num_sides
Geometric constants for all Tets.
Definition: cell_tet.h:74
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
Diagonal _diagonal_selection
The currently-selected diagonal used during refinement.
Definition: cell_tet.h:249
The libMesh namespace provides an interface to certain functionality in the library.
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: cell_tet10.C:87
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: cell_tet10.C:233
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_tet10.C:734
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: cell_tet10.C:77
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_tet10.h:265
virtual unsigned int n_sub_elem() const override
Definition: cell_tet10.h:95
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_vertices() const override final
Definition: cell_tet.h:86
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
static const unsigned short int _second_order_vertex_child_number[10]
Vector that names a child sharing each second order node.
Definition: cell_tet10.h:282
ElemMappingType mapping_type() const
Definition: elem.h:3120
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override
Definition: cell_tet10.C:113
virtual Real volume() const override
A specialization for computing the volume of a Tet10.
Definition: cell_tet10.C:622
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Builds a TRI6 built coincident with face i.
Definition: cell_tet10.C:204
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: cell_tet10.C:100
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
virtual unsigned int n_edges() const override final
Definition: cell_tet.h:91
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_tet10.C:564
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Definition: cell_tet10.C:182
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_tet10.C:776
static constexpr Real affine_tol
Default tolerance to use in has_affine_map().
Definition: elem.h:2061
static const int nodes_per_edge
Definition: cell_tet10.h:218
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 std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
Definition: cell_tet10.C:94
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
ElemType side_type(const unsigned int s) const override final
Definition: cell_tet10.C:790
virtual bool is_vertex(const unsigned int i) const override
Definition: cell_tet10.C:58
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_tet10.h:224
virtual Real volume() const
Definition: elem.C:3429
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 bool is_face(const unsigned int i) const override
Definition: cell_tet10.C:72
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
Definition: cell_tet10.C:401
static const unsigned short int _second_order_vertex_child_index[10]
Vector that names the child vertex index for each second order node.
Definition: cell_tet10.h:287
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
Definition: cell_tet10.C:193
virtual bool is_edge(const unsigned int i) const override
Definition: cell_tet10.C:65
static const int num_nodes
Geometric constants for Tet10.
Definition: cell_tet10.h:216
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
static const int nodes_per_side
Definition: cell_tet10.h:217