Line data Source code
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 :
45 : const unsigned int Tet14::side_nodes_map[Tet14::num_sides][Tet14::nodes_per_side] =
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 :
53 : const unsigned int Tet14::edge_nodes_map[Tet14::num_edges][Tet14::nodes_per_edge] =
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 64192497 : bool Tet14::is_vertex(const unsigned int i) const
67 : {
68 64192497 : if (i < 4)
69 18113662 : return true;
70 3399014 : return false;
71 : }
72 :
73 138240 : bool Tet14::is_edge(const unsigned int i) const
74 : {
75 138240 : if (i < 4 || i > 9)
76 55296 : return false;
77 6912 : return true;
78 : }
79 :
80 8291478 : bool Tet14::is_face(const unsigned int i) const
81 : {
82 8291478 : if (i > 9)
83 1610178 : return true;
84 545530 : return false;
85 : }
86 :
87 4019014 : bool Tet14::is_node_on_side(const unsigned int n,
88 : const unsigned int s) const
89 : {
90 318493 : libmesh_assert_less (s, n_sides());
91 318493 : return std::find(std::begin(side_nodes_map[s]),
92 4019014 : std::end(side_nodes_map[s]),
93 4019014 : n) != std::end(side_nodes_map[s]);
94 : }
95 :
96 : std::vector<unsigned>
97 92659868 : Tet14::nodes_on_side(const unsigned int s) const
98 : {
99 7721768 : libmesh_assert_less(s, n_sides());
100 100381636 : return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
101 : }
102 :
103 : std::vector<unsigned>
104 83370 : Tet14::nodes_on_edge(const unsigned int e) const
105 : {
106 6924 : libmesh_assert_less(e, n_edges());
107 159816 : return {std::begin(edge_nodes_map[e]), std::end(edge_nodes_map[e])};
108 : }
109 :
110 404268 : bool Tet14::is_node_on_edge(const unsigned int n,
111 : const unsigned int e) const
112 : {
113 27060 : libmesh_assert_less (e, n_edges());
114 27060 : return std::find(std::begin(edge_nodes_map[e]),
115 404268 : std::end(edge_nodes_map[e]),
116 404268 : 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 1611764 : 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 1611764 : 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 2064308 : 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 :
153 16690606 : bool Tet14::has_affine_map() const
154 : {
155 : // Make sure edges are straight
156 2775490 : Point v10 = this->point(1) - this->point(0);
157 18078351 : if (!v10.relative_fuzzy_equals
158 16690606 : ((this->point(4) - this->point(0))*2))
159 0 : return false;
160 2775490 : Point v21 = this->point(2) - this->point(1);
161 18078351 : if (!v21.relative_fuzzy_equals
162 16690606 : ((this->point(5) - this->point(1))*2))
163 0 : return false;
164 2775490 : Point v20 = this->point(2) - this->point(0);
165 18078351 : if (!v20.relative_fuzzy_equals
166 16690606 : ((this->point(6) - this->point(0))*2))
167 0 : return false;
168 2775490 : Point v30 = this->point(3) - this->point(0);
169 18078351 : if (!v30.relative_fuzzy_equals
170 16690606 : ((this->point(7) - this->point(0))*2))
171 0 : return false;
172 2775490 : Point v31 = this->point(3) - this->point(1);
173 18078351 : if (!v31.relative_fuzzy_equals
174 16690606 : ((this->point(8) - this->point(1))*2))
175 0 : return false;
176 2775490 : Point v32 = this->point(3) - this->point(2);
177 18078351 : if (!v32.relative_fuzzy_equals
178 16690606 : ((this->point(9) - this->point(2))*2))
179 0 : return false;
180 :
181 : // Make sure midface nodes are midface
182 16690606 : if (!(v20+v10).relative_fuzzy_equals
183 16690606 : ((this->point(10) - this->point(0))*3))
184 0 : return false;
185 16690606 : if (!(v30+v10).relative_fuzzy_equals
186 16690606 : ((this->point(11) - this->point(0))*3))
187 0 : return false;
188 16690606 : if (!(v31+v21).relative_fuzzy_equals
189 16690606 : ((this->point(12) - this->point(1))*3))
190 0 : return false;
191 16690606 : if (!(v30+v20).relative_fuzzy_equals
192 16690606 : ((this->point(13) - this->point(0))*3))
193 0 : return false;
194 :
195 1387745 : return true;
196 : }
197 :
198 :
199 :
200 118190377 : Order Tet14::default_order() const
201 : {
202 118190377 : return THIRD;
203 : }
204 :
205 :
206 :
207 82205214 : unsigned int Tet14::local_side_node(unsigned int side,
208 : unsigned int side_node) const
209 : {
210 6843676 : libmesh_assert_less (side, this->n_sides());
211 6843676 : libmesh_assert_less (side_node, Tet14::nodes_per_side);
212 :
213 82205214 : return Tet14::side_nodes_map[side][side_node];
214 : }
215 :
216 :
217 :
218 38869380 : unsigned int Tet14::local_edge_node(unsigned int edge,
219 : unsigned int edge_node) const
220 : {
221 2374854 : libmesh_assert_less (edge, this->n_edges());
222 2374854 : libmesh_assert_less (edge_node, Tet14::nodes_per_edge);
223 :
224 38869380 : return Tet14::edge_nodes_map[edge][edge_node];
225 : }
226 :
227 :
228 :
229 18844870 : std::unique_ptr<Elem> Tet14::build_side_ptr (const unsigned int i)
230 : {
231 18844870 : return this->simple_build_side_ptr<Tri7, Tet14>(i);
232 : }
233 :
234 :
235 :
236 13244036 : void Tet14::build_side_ptr (std::unique_ptr<Elem> & side,
237 : const unsigned int i)
238 : {
239 13244036 : this->simple_build_side_ptr<Tet14>(side, i, TRI7);
240 13244036 : }
241 :
242 :
243 :
244 5234640 : std::unique_ptr<Elem> Tet14::build_edge_ptr (const unsigned int i)
245 : {
246 5234640 : return this->simple_build_edge_ptr<Edge3,Tet14>(i);
247 : }
248 :
249 :
250 :
251 0 : void Tet14::build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i)
252 : {
253 0 : this->simple_build_edge_ptr<Tet14>(edge, i, EDGE3);
254 0 : }
255 :
256 :
257 :
258 0 : void Tet14::connectivity(const unsigned int sc,
259 : const IOPackage iop,
260 : std::vector<dof_id_type> & conn) const
261 : {
262 0 : libmesh_assert(_nodes);
263 0 : libmesh_assert_less (sc, this->n_sub_elem());
264 0 : libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
265 :
266 0 : switch (iop)
267 : {
268 0 : case TECPLOT:
269 : {
270 0 : conn.resize(8);
271 0 : switch (sc)
272 : {
273 :
274 :
275 : // Linear sub-tet 0
276 0 : case 0:
277 :
278 0 : conn[0] = this->node_id(0)+1;
279 0 : conn[1] = this->node_id(4)+1;
280 0 : conn[2] = this->node_id(6)+1;
281 0 : conn[3] = this->node_id(6)+1;
282 0 : conn[4] = this->node_id(7)+1;
283 0 : conn[5] = this->node_id(7)+1;
284 0 : conn[6] = this->node_id(7)+1;
285 0 : conn[7] = this->node_id(7)+1;
286 :
287 0 : return;
288 :
289 : // Linear sub-tet 1
290 0 : case 1:
291 :
292 0 : conn[0] = this->node_id(4)+1;
293 0 : conn[1] = this->node_id(1)+1;
294 0 : conn[2] = this->node_id(5)+1;
295 0 : conn[3] = this->node_id(5)+1;
296 0 : conn[4] = this->node_id(8)+1;
297 0 : conn[5] = this->node_id(8)+1;
298 0 : conn[6] = this->node_id(8)+1;
299 0 : conn[7] = this->node_id(8)+1;
300 :
301 0 : return;
302 :
303 : // Linear sub-tet 2
304 0 : case 2:
305 :
306 0 : conn[0] = this->node_id(5)+1;
307 0 : conn[1] = this->node_id(2)+1;
308 0 : conn[2] = this->node_id(6)+1;
309 0 : conn[3] = this->node_id(6)+1;
310 0 : conn[4] = this->node_id(9)+1;
311 0 : conn[5] = this->node_id(9)+1;
312 0 : conn[6] = this->node_id(9)+1;
313 0 : conn[7] = this->node_id(9)+1;
314 :
315 0 : return;
316 :
317 : // Linear sub-tet 3
318 0 : case 3:
319 :
320 0 : conn[0] = this->node_id(7)+1;
321 0 : conn[1] = this->node_id(8)+1;
322 0 : conn[2] = this->node_id(9)+1;
323 0 : conn[3] = this->node_id(9)+1;
324 0 : conn[4] = this->node_id(3)+1;
325 0 : conn[5] = this->node_id(3)+1;
326 0 : conn[6] = this->node_id(3)+1;
327 0 : conn[7] = this->node_id(3)+1;
328 :
329 0 : return;
330 :
331 : // Linear sub-tet 4
332 0 : case 4:
333 :
334 0 : conn[0] = this->node_id(4)+1;
335 0 : conn[1] = this->node_id(8)+1;
336 0 : conn[2] = this->node_id(6)+1;
337 0 : conn[3] = this->node_id(6)+1;
338 0 : conn[4] = this->node_id(7)+1;
339 0 : conn[5] = this->node_id(7)+1;
340 0 : conn[6] = this->node_id(7)+1;
341 0 : conn[7] = this->node_id(7)+1;
342 :
343 0 : return;
344 :
345 : // Linear sub-tet 5
346 0 : case 5:
347 :
348 0 : conn[0] = this->node_id(4)+1;
349 0 : conn[1] = this->node_id(5)+1;
350 0 : conn[2] = this->node_id(6)+1;
351 0 : conn[3] = this->node_id(6)+1;
352 0 : conn[4] = this->node_id(8)+1;
353 0 : conn[5] = this->node_id(8)+1;
354 0 : conn[6] = this->node_id(8)+1;
355 0 : conn[7] = this->node_id(8)+1;
356 :
357 0 : return;
358 :
359 : // Linear sub-tet 6
360 0 : case 6:
361 :
362 0 : conn[0] = this->node_id(5)+1;
363 0 : conn[1] = this->node_id(9)+1;
364 0 : conn[2] = this->node_id(6)+1;
365 0 : conn[3] = this->node_id(6)+1;
366 0 : conn[4] = this->node_id(8)+1;
367 0 : conn[5] = this->node_id(8)+1;
368 0 : conn[6] = this->node_id(8)+1;
369 0 : conn[7] = this->node_id(8)+1;
370 :
371 0 : return;
372 :
373 : // Linear sub-tet 7
374 0 : case 7:
375 :
376 0 : conn[0] = this->node_id(7)+1;
377 0 : conn[1] = this->node_id(6)+1;
378 0 : conn[2] = this->node_id(9)+1;
379 0 : conn[3] = this->node_id(9)+1;
380 0 : conn[4] = this->node_id(8)+1;
381 0 : conn[5] = this->node_id(8)+1;
382 0 : conn[6] = this->node_id(8)+1;
383 0 : conn[7] = this->node_id(8)+1;
384 :
385 0 : return;
386 :
387 :
388 0 : default:
389 0 : libmesh_error_msg("Invalid sc = " << sc);
390 : }
391 : }
392 :
393 0 : 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 0 : conn.resize(Tet14::num_nodes);
400 0 : for (auto i : index_range(conn))
401 0 : conn[i] = this->node_id(i);
402 0 : return;
403 : }
404 :
405 0 : default:
406 0 : libmesh_error_msg("Unsupported IO package " << iop);
407 : }
408 : }
409 :
410 :
411 :
412 23804880 : unsigned int Tet14::n_second_order_adjacent_vertices (const unsigned int n) const
413 : {
414 23804880 : switch (n)
415 : {
416 402336 : case 4:
417 : case 5:
418 : case 6:
419 : case 7:
420 : case 8:
421 : case 9:
422 402336 : return 2;
423 :
424 9521952 : case 10:
425 : case 11:
426 : case 12:
427 : case 13:
428 9521952 : return 3;
429 :
430 0 : default:
431 0 : 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 0 : Tet14::second_order_child_vertex (const unsigned int n) const
458 : {
459 0 : libmesh_assert_greater_equal (n, this->n_vertices());
460 0 : libmesh_assert_less (n, this->n_nodes());
461 0 : return std::pair<unsigned short int, unsigned short int>
462 0 : (_second_order_vertex_child_number[n],
463 0 : _second_order_vertex_child_index[n]);
464 : }
465 :
466 :
467 :
468 57131712 : unsigned short int Tet14::second_order_adjacent_vertex (const unsigned int n,
469 : const unsigned int v) const
470 : {
471 1609344 : libmesh_assert_greater_equal (n, this->n_vertices());
472 1609344 : libmesh_assert_less (n, this->n_nodes());
473 1609344 : libmesh_assert_less (v, 3);
474 1609344 : libmesh_assert (n > 9 || v < 2); // Only face nodes have multiple adjacencies
475 57131712 : 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 :
500 : const Real Tet14::_embedding_matrix[Tet14::num_children][Tet14::num_nodes][Tet14::num_nodes] =
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 33726609 : 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 33726609 : this->choose_diagonal();
663 :
664 : // Permuted j and k indices
665 : unsigned int
666 2039726 : jp=j,
667 2039726 : kp=k;
668 :
669 33726609 : if ((i>3) && (this->_diagonal_selection!=DIAG_02_13))
670 : {
671 : // Just the enum value cast to an unsigned int...
672 246236 : 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 2474852 : 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 2474852 : 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 2474852 : 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 33726609 : return this->_embedding_matrix[i][jp][kp];
711 : }
712 :
713 :
714 : const std::vector<std::pair<unsigned char, unsigned char>> &
715 7634568 : Tet14::parent_bracketing_nodes(unsigned int c,
716 : unsigned int n) const
717 : {
718 : // Choose an optimal diagonal, if one has not already been selected
719 7634568 : this->choose_diagonal();
720 :
721 : // Just the enum value cast to an unsigned int...
722 7634568 : const unsigned ds = static_cast<unsigned int>(this->_diagonal_selection);
723 :
724 7634568 : return Tet14::_parent_bracketing_nodes[ds][c][n];
725 : }
726 :
727 :
728 : const std::vector<std::pair<unsigned char, unsigned char>>
729 : Tet14::_parent_bracketing_nodes[3][Tet14::num_children][Tet14::num_nodes] =
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 950605 : void Tet14::permute(unsigned int perm_num)
813 : {
814 60480 : libmesh_assert_less (perm_num, 12);
815 :
816 950605 : const unsigned int side = perm_num % 4;
817 950605 : const unsigned int rotate = perm_num / 4;
818 :
819 1912888 : for (unsigned int i = 0; i != rotate; ++i)
820 : {
821 962283 : swap3nodes(0,1,2);
822 901607 : swap3nodes(4,5,6);
823 901607 : swap3nodes(7,8,9);
824 901607 : swap3nodes(11,12,13);
825 901607 : swap3neighbors(1,2,3);
826 : }
827 :
828 950605 : switch (side) {
829 15898 : case 0:
830 15898 : break;
831 221957 : case 1:
832 221957 : swap3nodes(0,2,3);
833 207597 : swap3nodes(4,5,8);
834 207597 : swap3nodes(6,9,7);
835 207597 : swap3nodes(10,12,11);
836 207597 : swap3neighbors(0,2,1);
837 207597 : break;
838 231471 : case 2:
839 231471 : swap3nodes(2,0,3);
840 216415 : swap3nodes(5,4,8);
841 216415 : swap3nodes(6,7,9);
842 216415 : swap3nodes(10,11,12);
843 216415 : swap3neighbors(0,1,2);
844 216415 : break;
845 240354 : case 3:
846 240354 : swap3nodes(2,1,3);
847 225188 : swap3nodes(5,8,9);
848 225188 : swap3nodes(6,4,7);
849 225188 : swap3nodes(10,11,13);
850 225188 : swap3neighbors(0,1,3);
851 225188 : break;
852 0 : default:
853 0 : libmesh_error();
854 : }
855 950605 : }
856 :
857 :
858 6912 : void Tet14::flip(BoundaryInfo * boundary_info)
859 : {
860 576 : libmesh_assert(boundary_info);
861 :
862 6912 : swap2nodes(0,2);
863 6912 : swap2nodes(4,5);
864 6912 : swap2nodes(7,9);
865 6912 : swap2nodes(11,12);
866 576 : swap2neighbors(1,2);
867 6912 : swap2boundarysides(1,2,boundary_info);
868 6912 : swap2boundaryedges(0,1,boundary_info);
869 6912 : swap2boundaryedges(3,5,boundary_info);
870 6912 : }
871 :
872 :
873 27648 : ElemType Tet14::side_type (const unsigned int libmesh_dbg_var(s)) const
874 : {
875 2304 : libmesh_assert_less (s, 4);
876 27648 : return TRI7;
877 : }
878 :
879 :
880 : } // namespace libMesh
|