libMesh
cell_hex8.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_hex8.h"
21 #include "libmesh/edge_edge2.h"
22 #include "libmesh/face_quad4.h"
23 #include "libmesh/enum_io_package.h"
24 #include "libmesh/enum_order.h"
25 #include "libmesh/fe_lagrange_shape_1D.h"
26 
27 // C++ includes
28 #include <array>
29 
30 namespace libMesh
31 {
32 
33 
34 
35 
36 // ------------------------------------------------------------
37 // Hex8 class static member initializations
38 const int Hex8::num_nodes;
39 const int Hex8::nodes_per_side;
40 const int Hex8::nodes_per_edge;
41 
43  {
44  {0, 3, 2, 1}, // Side 0
45  {0, 1, 5, 4}, // Side 1
46  {1, 2, 6, 5}, // Side 2
47  {2, 3, 7, 6}, // Side 3
48  {3, 0, 4, 7}, // Side 4
49  {4, 5, 6, 7} // Side 5
50  };
51 
53  {
54  {0, 1}, // Edge 0
55  {1, 2}, // Edge 1
56  {2, 3}, // Edge 2
57  {0, 3}, // Edge 3
58  {0, 4}, // Edge 4
59  {1, 5}, // Edge 5
60  {2, 6}, // Edge 6
61  {3, 7}, // Edge 7
62  {4, 5}, // Edge 8
63  {5, 6}, // Edge 9
64  {6, 7}, // Edge 10
65  {4, 7} // Edge 11
66  };
67 
68 // ------------------------------------------------------------
69 // Hex8 class member functions
70 
71 bool Hex8::is_vertex(const unsigned int libmesh_dbg_var(n)) const
72 {
73  libmesh_assert_not_equal_to (n, invalid_uint);
74  return true;
75 }
76 
77 bool Hex8::is_edge(const unsigned int) const
78 {
79  return false;
80 }
81 
82 bool Hex8::is_face(const unsigned int) const
83 {
84  return false;
85 }
86 
87 bool Hex8::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 Hex8::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 Hex8::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 Hex8::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 
122 {
123  // Make sure x-edge endpoints are affine
124  Point v = this->point(1) - this->point(0);
125  if (!v.relative_fuzzy_equals(this->point(2) - this->point(3), affine_tol) ||
126  !v.relative_fuzzy_equals(this->point(5) - this->point(4), affine_tol) ||
127  !v.relative_fuzzy_equals(this->point(6) - this->point(7), affine_tol))
128  return false;
129  // Make sure xz-faces are identical parallelograms
130  v = this->point(4) - this->point(0);
131  if (!v.relative_fuzzy_equals(this->point(7) - this->point(3), affine_tol))
132  return false;
133  // If all the above checks out, the map is affine
134  return true;
135 }
136 
137 
138 
140 {
141  return FIRST;
142 }
143 
144 
145 
146 std::unique_ptr<Elem> Hex8::build_side_ptr (const unsigned int i)
147 {
148  return this->simple_build_side_ptr<Quad4, Hex8>(i);
149 }
150 
151 
152 
153 void Hex8::build_side_ptr (std::unique_ptr<Elem> & side,
154  const unsigned int i)
155 {
156  this->simple_build_side_ptr<Hex8>(side, i, QUAD4);
157 }
158 
159 
160 
161 std::unique_ptr<Elem> Hex8::build_edge_ptr (const unsigned int i)
162 {
163  return this->simple_build_edge_ptr<Edge2,Hex8>(i);
164 }
165 
166 
167 
168 void Hex8::build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i)
169 {
170  this->simple_build_edge_ptr<Hex8>(edge, i, EDGE2);
171 }
172 
173 
174 
175 void Hex8::connectivity(const unsigned int libmesh_dbg_var(sc),
176  const IOPackage iop,
177  std::vector<dof_id_type> & conn) const
178 {
180  libmesh_assert_less (sc, this->n_sub_elem());
181  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
182 
183  conn.resize(8);
184 
185  switch (iop)
186  {
187  case TECPLOT:
188  {
189  conn[0] = this->node_id(0)+1;
190  conn[1] = this->node_id(1)+1;
191  conn[2] = this->node_id(2)+1;
192  conn[3] = this->node_id(3)+1;
193  conn[4] = this->node_id(4)+1;
194  conn[5] = this->node_id(5)+1;
195  conn[6] = this->node_id(6)+1;
196  conn[7] = this->node_id(7)+1;
197  return;
198  }
199 
200  case VTK:
201  {
202  for (auto i : index_range(conn))
203  conn[i] = this->node_id(i);
204  return;
205  }
206 
207  default:
208  libmesh_error_msg("Unsupported IO package " << iop);
209  }
210 }
211 
212 
213 
214 #ifdef LIBMESH_ENABLE_AMR
215 
217  {
218  // The 8 children of the Hex-type elements can be thought of as being
219  // associated with the 8 vertices of the Hex. Some of the children are
220  // numbered the same as their corresponding vertex, while some are
221  // not. The children which are numbered differently have been marked
222  // with ** in the comments below.
223 
224  // embedding matrix for child 0 (child 0 is associated with vertex 0)
225  {
226  // 0 1 2 3 4 5 6 7
227  { 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0
228  { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 1
229  { .25, .25, .25, .25, 0.0, 0.0, 0.0, 0.0}, // 2
230  { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0}, // 3
231  { 0.5, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0}, // 4
232  { .25, .25, 0.0, 0.0, .25, .25, 0.0, 0.0}, // 5
233  {.125, .125, .125, .125, .125, .125, .125, .125}, // 6
234  { .25, 0.0, 0.0, .25, .25, 0.0, 0.0, .25} // 7
235  },
236 
237  // embedding matrix for child 1 (child 1 is associated with vertex 1)
238  {
239  // 0 1 2 3 4 5 6 7
240  { 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 0
241  { 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 1
242  { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0}, // 2
243  { .25, .25, .25, .25, 0.0, 0.0, 0.0, 0.0}, // 3
244  { .25, .25, 0.0, 0.0, .25, .25, 0.0, 0.0}, // 4
245  { 0.0, 0.5, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0}, // 5
246  { 0.0, .25, .25, 0.0, 0.0, .25, .25, 0.0}, // 6
247  {.125, .125, .125, .125, .125, .125, .125, .125} // 7
248  },
249 
250  // embedding matrix for child 2 (child 2 is associated with vertex 3**)
251  {
252  // 0 1 2 3 4 5 6 7
253  { 0.5, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0}, // 0
254  { .25, .25, .25, .25, 0.0, 0.0, 0.0, 0.0}, // 1
255  { 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 2
256  { 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0}, // 3
257  { .25, 0.0, 0.0, .25, .25, 0.0, 0.0, .25}, // 4
258  {.125, .125, .125, .125, .125, .125, .125, .125}, // 5
259  { 0.0, 0.0, .25, .25, 0.0, 0.0, .25, .25}, // 6
260  { 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5} // 7
261  },
262 
263  // embedding matrix for child 3 (child 3 is associated with vertex 2**)
264  {
265  // 0 1 2 3 4 5 6 7
266  { .25, .25, .25, .25, 0.0, 0.0, 0.0, 0.0}, // 0
267  { 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0}, // 1
268  { 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0}, // 2
269  { 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0}, // 3
270  {.125, .125, .125, .125, .125, .125, .125, .125}, // 4
271  { 0.0, .25, .25, 0.0, 0.0, .25, .25, 0.0}, // 5
272  { 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5, 0.0}, // 6
273  { 0.0, 0.0, .25, .25, 0.0, 0.0, .25, .25} // 7
274  },
275 
276  // embedding matrix for child 4 (child 4 is associated with vertex 4)
277  {
278  // 0 1 2 3 4 5 6 7
279  { 0.5, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0}, // 0
280  { .25, .25, 0.0, 0.0, .25, .25, 0.0, 0.0}, // 1
281  {.125, .125, .125, .125, .125, .125, .125, .125}, // 2
282  { .25, 0.0, 0.0, .25, .25, 0.0, 0.0, .25}, // 3
283  { 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0}, // 4
284  { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0}, // 5
285  { 0.0, 0.0, 0.0, 0.0, .25, .25, .25, .25}, // 6
286  { 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.5} // 7
287  },
288 
289  // embedding matrix for child 5 (child 5 is associated with vertex 5)
290  {
291  // 0 1 2 3 4 5 6 7
292  { .25, .25, 0.0, 0.0, .25, .25, 0.0, 0.0}, // 0
293  { 0.0, 0.5, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0}, // 1
294  { 0.0, .25, .25, 0.0, 0.0, .25, .25, 0.0}, // 2
295  {.125, .125, .125, .125, .125, .125, .125, .125}, // 3
296  { 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0}, // 4
297  { 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0}, // 5
298  { 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 6
299  { 0.0, 0.0, 0.0, 0.0, .25, .25, .25, .25} // 7
300  },
301 
302  // embedding matrix for child 6 (child 6 is associated with vertex 7**)
303  {
304  // 0 1 2 3 4 5 6 7
305  { .25, 0.0, 0.0, .25, .25, 0.0, 0.0, .25}, // 0
306  {.125, .125, .125, .125, .125, .125, .125, .125}, // 1
307  { 0.0, 0.0, .25, .25, 0.0, 0.0, .25, .25}, // 2
308  { 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5}, // 3
309  { 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.5}, // 4
310  { 0.0, 0.0, 0.0, 0.0, .25, .25, .25, .25}, // 5
311  { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5}, // 6
312  { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0} // 7
313  },
314 
315  // embedding matrix for child 7 (child 7 is associated with vertex 6**)
316  {
317  // 0 1 2 3 4 5 6 7
318  {.125, .125, .125, .125, .125, .125, .125, .125}, // 0
319  { 0.0, .25, .25, 0.0, 0.0, .25, .25, 0.0}, // 1
320  { 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.5, 0.0}, // 2
321  { 0.0, 0.0, .25, .25, 0.0, 0.0, .25, .25}, // 3
322  { 0.0, 0.0, 0.0, 0.0, .25, .25, .25, .25}, // 4
323  { 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0}, // 5
324  { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0}, // 6
325  { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5} // 7
326  }
327  };
328 
329 
330 
331 
332 #endif
333 
334 
335 
337  const Point & x0, const Point & x1, const Point & x2, const Point & x3,
338  const Point & x4, const Point & x5, const Point & x6, const Point & x7)
339 {
340  // The Jacobian is dx/d(xi) dot (dx/d(eta) cross dx/d(zeta)), where
341  // dx/d(xi) = a1*eta*zeta + b1*eta + c1*zeta + d1
342  // dx/d(eta) = a2*xi*zeta + b2*xi + c2*zeta + d2
343  // dx/d(zeta) = a3*xi*eta + b3*xi + c3*eta + d3
344 
345  // Notes:
346  // 1.) Several of these coefficient vectors are equal, as noted below.
347  // 2.) These are all off by a factor of 8, but this cancels when we
348  // divide by the volume, which will also be off by the same
349  // factor.
350  Point
351  a1 = -x0 + x1 - x2 + x3 + x4 - x5 + x6 - x7,
352  a2 = a1,
353  a3 = a1;
354 
355  Point
356  b1 = x0 - x1 + x2 - x3 + x4 - x5 + x6 - x7,
357  b2 = b1,
358  b3 = x0 - x1 - x2 + x3 - x4 + x5 + x6 - x7;
359 
360  Point
361  c1 = b3,
362  c2 = x0 + x1 - x2 - x3 - x4 - x5 + x6 + x7,
363  c3 = c2;
364 
365  Point
366  d1 = -x0 + x1 + x2 - x3 - x4 + x5 + x6 - x7,
367  d2 = -x0 - x1 + x2 + x3 - x4 - x5 + x6 + x7,
368  d3 = -x0 - x1 - x2 - x3 + x4 + x5 + x6 + x7;
369 
370  // Use 2x2x2 quadrature to compute the integral of each basis
371  // function (as defined on the [-1,1]^3 reference domain). We use
372  // a quadrature rule which is exact for tri-cubics. The weights for
373  // this rule are all equal to 1.
374  static const Real q[2] = {-std::sqrt(Real(3))/3, std::sqrt(Real(3))/3.};
375 
376  // Indices for computing tensor product basis functions. This is
377  // copied from fe_lagrange_shape_3D.C
378  static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0};
379  static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1};
380  static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1};
381 
382  // Compute nodal volumes
383  std::array<Real, Hex8::num_nodes> V{};
384 
385  for (const auto & xi : q)
386  for (const auto & eta : q)
387  for (const auto & zeta : q)
388  {
389  Real jxw = triple_product(a1*eta*zeta + b1*eta + c1*zeta + d1,
390  a2*xi*zeta + b2*xi + c2*zeta + d2,
391  a3*xi*eta + b3*xi + c3*eta + d3);
392 
393  for (int i=0; i<Hex8::num_nodes; ++i)
394  V[i] += jxw *
395  fe_lagrange_1D_linear_shape(i0[i], xi) *
396  fe_lagrange_1D_linear_shape(i1[i], eta) *
397  fe_lagrange_1D_linear_shape(i2[i], zeta);
398  }
399 
400  // Compute centroid
401  return
402  (x0*V[0] + x1*V[1] + x2*V[2] + x3*V[3] + x4*V[4] + x5*V[5] + x6*V[6] + x7*V[7]) /
403  (V[0] + V[1] + V[2] + V[3] + V[4] + V[5] + V[6] + V[7]);
404 }
405 
406 
408 {
410  (point(0), point(1), point(2), point(3),
411  point(4), point(5), point(6), point(7));
412 }
413 
414 
416 {
417  // Make copies of our points. It makes the subsequent calculations a bit
418  // shorter and avoids dereferencing the same pointer multiple times.
419  Point
420  x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3),
421  x4 = point(4), x5 = point(5), x6 = point(6), x7 = point(7);
422 
423  // Construct constant data vectors. The notation is:
424  // \vec{x}_{\xi} = \vec{a1}*eta*zeta + \vec{b1}*eta + \vec{c1}*zeta + \vec{d1}
425  // \vec{x}_{\eta} = \vec{a2}*xi*zeta + \vec{b2}*xi + \vec{c2}*zeta + \vec{d2}
426  // \vec{x}_{\zeta} = \vec{a3}*xi*eta + \vec{b3}*xi + \vec{c3}*eta + \vec{d3}
427  // but it turns out that a1, a2, and a3 are not needed for the volume calculation.
428 
429  // Build up the 6 unique vectors which make up dx/dxi, dx/deta, and dx/dzeta.
430  Point q[6] =
431  {
432  /*b1*/ x0 - x1 + x2 - x3 + x4 - x5 + x6 - x7, /*=b2*/
433  /*c1*/ x0 - x1 - x2 + x3 - x4 + x5 + x6 - x7, /*=b3*/
434  /*d1*/ -x0 + x1 + x2 - x3 - x4 + x5 + x6 - x7,
435  /*c2*/ x0 + x1 - x2 - x3 - x4 - x5 + x6 + x7, /*=c3*/
436  /*d2*/ -x0 - x1 + x2 + x3 - x4 - x5 + x6 + x7,
437  /*d3*/ -x0 - x1 - x2 - x3 + x4 + x5 + x6 + x7
438  };
439 
440  // We could check for a linear element, but it's probably faster to
441  // just compute the result...
442  return
443  (triple_product(q[0], q[4], q[3]) +
444  triple_product(q[2], q[0], q[1]) +
445  triple_product(q[1], q[3], q[5])) / 192. +
446  triple_product(q[2], q[4], q[5]) / 64.;
447 }
448 
451 {
452  return Elem::loose_bounding_box();
453 }
454 
455 
456 void
457 Hex8::permute(unsigned int perm_num)
458 {
459  libmesh_assert_less (perm_num, 24);
460  const unsigned int side = perm_num % 6;
461  const unsigned int rotate = perm_num / 6;
462 
463  for (unsigned int i = 0; i != rotate; ++i)
464  {
465  swap4nodes(0,1,2,3);
466  swap4nodes(4,5,6,7);
467  swap4neighbors(1,2,3,4);
468  }
469 
470  switch (side) {
471  case 0:
472  break;
473  case 1:
474  swap4nodes(3,7,4,0);
475  swap4nodes(2,6,5,1);
476  swap4neighbors(0,3,5,1);
477  break;
478  case 2:
479  swap4nodes(0,4,5,1);
480  swap4nodes(3,7,6,2);
481  swap4neighbors(0,4,5,2);
482  break;
483  case 3:
484  swap4nodes(0,4,7,3);
485  swap4nodes(1,5,6,2);
486  swap4neighbors(0,1,5,3);
487  break;
488  case 4:
489  swap4nodes(1,5,4,0);
490  swap4nodes(2,6,7,3);
491  swap4neighbors(0,2,5,4);
492  break;
493  case 5:
494  swap2nodes(0,7);
495  swap2nodes(1,6);
496  swap2nodes(2,5);
497  swap2nodes(3,4);
498  swap2neighbors(0,5);
499  swap2neighbors(1,3);
500  break;
501  default:
502  libmesh_error();
503  }
504 }
505 
506 
507 void
508 Hex8::flip(BoundaryInfo * boundary_info)
509 {
510  libmesh_assert(boundary_info);
511 
512  swap2nodes(0,1);
513  swap2nodes(2,3);
514  swap2nodes(4,5);
515  swap2nodes(6,7);
516  swap2neighbors(2,4);
517  swap2boundarysides(2,4,boundary_info);
518  swap2boundaryedges(1,3,boundary_info);
519  swap2boundaryedges(4,5,boundary_info);
520  swap2boundaryedges(6,7,boundary_info);
521  swap2boundaryedges(9,11,boundary_info);
522 }
523 
524 
525 ElemType
526 Hex8::side_type (const unsigned int libmesh_dbg_var(s)) const
527 {
528  libmesh_assert_less (s, 6);
529  return QUAD4;
530 }
531 
532 
533 } // namespace libMesh
static const int num_nodes
Geometric constants for Hex8.
Definition: cell_hex8.h:165
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: cell_hex8.C:97
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE2 built coincident with edge i.
Definition: cell_hex8.C:161
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 std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Builds a QUAD4 built coincident with face i.
Definition: cell_hex8.C:146
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
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: cell_hex8.C:175
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
virtual BoundingBox loose_bounding_box() const
Definition: elem.C:3465
virtual Real volume() const override
A specialization for computing the area of a hexahedron with flat sides.
Definition: cell_hex8.C:415
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
Definition: cell_hex8.C:104
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
The libMesh namespace provides an interface to certain functionality in the library.
virtual unsigned int n_sub_elem() const override
Definition: cell_hex8.h:86
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
Definition: cell_hex8.h:173
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_hex8.h:240
ElemType side_type(const unsigned int s) const override final
Definition: cell_hex8.C:526
virtual Order default_order() const override
Definition: cell_hex8.C:139
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1029
void swap4nodes(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four node_ptrs, "rotating" them.
Definition: elem.h:2143
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
static const int num_edges
Definition: cell_hex.h:74
static const int num_children
Definition: cell_hex.h:75
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_hex8.h:179
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
static const int nodes_per_edge
Definition: cell_hex8.h:167
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_hex8.C:457
static constexpr Real affine_tol
Default tolerance to use in has_affine_map().
Definition: elem.h:2061
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: cell_hex8.C:87
virtual bool has_affine_map() const override
Definition: cell_hex8.C:121
RealTensorValue rotate(MeshBase &mesh, const Real phi, const Real theta=0., const Real psi=0.)
Rotates the mesh in the xy plane.
virtual unsigned int n_edges() const override final
Definition: cell_hex.h:90
virtual bool is_face(const unsigned int i) const override
Definition: cell_hex8.C:82
virtual Point true_centroid() const override
We compute the centroid of the Hex using a customized numerical quadrature approach that avoids unnec...
Definition: cell_hex8.C:407
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:2102
static const int nodes_per_side
Definition: cell_hex8.h:166
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
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_hex8.C:508
virtual bool is_edge(const unsigned int i) const override
Definition: cell_hex8.C:77
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
Definition: cell_hex8.C:110
virtual unsigned int n_sides() const override final
Definition: cell_hex.h:80
void swap4neighbors(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four neighbor_ptrs, "rotating" them.
Definition: elem.h:2153
virtual BoundingBox loose_bounding_box() const override
Builds a bounding box out of the nodal positions.
Definition: cell_hex8.C:450
static const int num_sides
Geometric constants for all Hexes.
Definition: cell_hex.h:73
static Point centroid_from_points(const Point &x0, const Point &x1, const Point &x2, const Point &x3, const Point &x4, const Point &x5, const Point &x6, const Point &x7)
Class static helper function that computes the centroid of a hexahedral region from a set of input po...
Definition: cell_hex8.C:336
virtual bool is_vertex(const unsigned int i) const override
Definition: cell_hex8.C:71
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
Real fe_lagrange_1D_linear_shape(const unsigned int i, const Real xi)