libMesh
face_quad8.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 // Local includes
19 #include "libmesh/edge_edge3.h"
20 #include "libmesh/face_quad8.h"
21 #include "libmesh/enum_io_package.h"
22 #include "libmesh/enum_order.h"
23 
24 namespace libMesh
25 {
26 
27 
28 
29 
30 // ------------------------------------------------------------
31 // Quad8 class static member initializations
32 const int Quad8::num_nodes;
33 const int Quad8::nodes_per_side;
34 
36  {
37  {0, 1, 4}, // Side 0
38  {1, 2, 5}, // Side 1
39  {2, 3, 6}, // Side 2
40  {3, 0, 7} // Side 3
41  };
42 
43 
44 #ifdef LIBMESH_ENABLE_AMR
45 
47  {
48  // embedding matrix for child 0
49  {
50  // 0 1 2 3 4 5 6 7
51  { 1.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 0
52  { 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000 }, // 1
53  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 2
54  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000 }, // 3
55  { 0.375000, -0.125000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000, 0.00000 }, // 4
56  { -0.187500, -0.187500, -0.187500, -0.187500, 0.750000, 0.375000, 0.250000, 0.375000 }, // 5
57  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.250000, 0.375000, 0.750000 }, // 6
58  { 0.375000, 0.00000, 0.00000, -0.125000, 0.00000, 0.00000, 0.00000, 0.750000 } // 7
59  },
60 
61  // embedding matrix for child 1
62  {
63  // 0 1 2 3 4 5 6 7
64  { 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000 }, // 0
65  { 0.00000, 1.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 1
66  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000 }, // 2
67  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 3
68  { -0.125000, 0.375000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000, 0.00000 }, // 4
69  { 0.00000, 0.375000, -0.125000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000 }, // 5
70  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.750000, 0.375000, 0.250000 }, // 6
71  { -0.187500, -0.187500, -0.187500, -0.187500, 0.750000, 0.375000, 0.250000, 0.375000 } // 7
72  },
73 
74  // embedding matrix for child 2
75  {
76  // 0 1 2 3 4 5 6 7
77  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000 }, // 0
78  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 1
79  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000 }, // 2
80  { 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 3
81  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.250000, 0.375000, 0.750000 }, // 4
82  { -0.187500, -0.187500, -0.187500, -0.187500, 0.250000, 0.375000, 0.750000, 0.375000 }, // 5
83  { 0.00000, 0.00000, -0.125000, 0.375000, 0.00000, 0.00000, 0.750000, 0.00000 }, // 6
84  { -0.125000, 0.00000, 0.00000, 0.375000, 0.00000, 0.00000, 0.00000, 0.750000 } // 7
85  },
86 
87  // embedding matrix for child 3
88  {
89  // 0 1 2 3 4 5 6 7
90  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 0
91  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000 }, // 1
92  { 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 2
93  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000 }, // 3
94  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.750000, 0.375000, 0.250000 }, // 4
95  { 0.00000, -0.125000, 0.375000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000 }, // 5
96  { 0.00000, 0.00000, 0.375000, -0.125000, 0.00000, 0.00000, 0.750000, 0.00000 }, // 6
97  { -0.187500, -0.187500, -0.187500, -0.187500, 0.250000, 0.375000, 0.750000, 0.375000 } // 7
98  }
99  };
100 
101 
102 #endif
103 
104 
105 // ------------------------------------------------------------
106 // Quad8 class member functions
107 
108 bool Quad8::is_vertex(const unsigned int i) const
109 {
110  if (i < 4)
111  return true;
112  return false;
113 }
114 
115 bool Quad8::is_edge(const unsigned int i) const
116 {
117  if (i < 4)
118  return false;
119  return true;
120 }
121 
122 bool Quad8::is_face(const unsigned int) const
123 {
124  return false;
125 }
126 
127 bool Quad8::is_node_on_side(const unsigned int n,
128  const unsigned int s) const
129 {
130  libmesh_assert_less (s, n_sides());
131  return std::find(std::begin(side_nodes_map[s]),
132  std::end(side_nodes_map[s]),
133  n) != std::end(side_nodes_map[s]);
134 }
135 
136 std::vector<unsigned>
137 Quad8::nodes_on_side(const unsigned int s) const
138 {
139  libmesh_assert_less(s, n_sides());
140  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
141 }
142 
143 std::vector<unsigned>
144 Quad8::nodes_on_edge(const unsigned int e) const
145 {
146  return nodes_on_side(e);
147 }
148 
150 {
151  // make sure corners form a parallelogram
152  Point v = this->point(1) - this->point(0);
153  if (!v.relative_fuzzy_equals(this->point(2) - this->point(3), affine_tol))
154  return false;
155  // make sure sides are straight
156  v /= 2;
157  if (!v.relative_fuzzy_equals(this->point(4) - this->point(0), affine_tol) ||
158  !v.relative_fuzzy_equals(this->point(6) - this->point(3), affine_tol))
159  return false;
160  v = (this->point(3) - this->point(0))/2;
161  if (!v.relative_fuzzy_equals(this->point(7) - this->point(0), affine_tol) ||
162  !v.relative_fuzzy_equals(this->point(5) - this->point(1), affine_tol))
163  return false;
164  return true;
165 }
166 
167 
168 
170 {
171  return SECOND;
172 }
173 
174 
175 
176 dof_id_type Quad8::key (const unsigned int s) const
177 {
178  libmesh_assert_less (s, this->n_sides());
179 
180  switch (s)
181  {
182  case 0:
183 
184  return
185  this->compute_key (this->node_id(4));
186 
187  case 1:
188 
189  return
190  this->compute_key (this->node_id(5));
191 
192  case 2:
193 
194  return
195  this->compute_key (this->node_id(6));
196 
197  case 3:
198 
199  return
200  this->compute_key (this->node_id(7));
201 
202  default:
203  libmesh_error_msg("Invalid side s = " << s);
204  }
205 }
206 
207 
208 
209 unsigned int Quad8::local_side_node(unsigned int side,
210  unsigned int side_node) const
211 {
212  libmesh_assert_less (side, this->n_sides());
213  libmesh_assert_less (side_node, Quad8::nodes_per_side);
214 
215  return Quad8::side_nodes_map[side][side_node];
216 }
217 
218 
219 
220 std::unique_ptr<Elem> Quad8::build_side_ptr (const unsigned int i)
221 {
222  return this->simple_build_side_ptr<Edge3, Quad8>(i);
223 }
224 
225 
226 
227 void Quad8::build_side_ptr (std::unique_ptr<Elem> & side,
228  const unsigned int i)
229 {
230  this->simple_build_side_ptr<Quad8>(side, i, EDGE3);
231 }
232 
233 
234 
235 
236 
237 
238 void Quad8::connectivity(const unsigned int sf,
239  const IOPackage iop,
240  std::vector<dof_id_type> & conn) const
241 {
242  libmesh_assert_less (sf, this->n_sub_elem());
243  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
244 
245  switch (iop)
246  {
247  // Note: TECPLOT connectivity is output as four triangles with
248  // a central quadrilateral. Therefore, the first four connectivity
249  // arrays are degenerate quads (triangles in Tecplot).
250  case TECPLOT:
251  {
252  // Create storage
253  conn.resize(4);
254 
255  switch(sf)
256  {
257  case 0:
258  // linear sub-tri 0
259  conn[0] = this->node_id(0)+1;
260  conn[1] = this->node_id(4)+1;
261  conn[2] = this->node_id(7)+1;
262  conn[3] = this->node_id(7)+1;
263 
264  return;
265 
266  case 1:
267  // linear sub-tri 1
268  conn[0] = this->node_id(4)+1;
269  conn[1] = this->node_id(1)+1;
270  conn[2] = this->node_id(5)+1;
271  conn[3] = this->node_id(5)+1;
272 
273  return;
274 
275  case 2:
276  // linear sub-tri 2
277  conn[0] = this->node_id(5)+1;
278  conn[1] = this->node_id(2)+1;
279  conn[2] = this->node_id(6)+1;
280  conn[3] = this->node_id(6)+1;
281 
282  return;
283 
284  case 3:
285  // linear sub-tri 3
286  conn[0] = this->node_id(7)+1;
287  conn[1] = this->node_id(6)+1;
288  conn[2] = this->node_id(3)+1;
289  conn[3] = this->node_id(3)+1;
290 
291  return;
292 
293  case 4:
294  // linear sub-quad
295  conn[0] = this->node_id(4)+1;
296  conn[1] = this->node_id(5)+1;
297  conn[2] = this->node_id(6)+1;
298  conn[3] = this->node_id(7)+1;
299 
300  return;
301 
302  default:
303  libmesh_error_msg("Invalid sf = " << sf);
304  }
305  }
306 
307 
308  // VTK connectivity for this element matches libmesh's own.
309  case VTK:
310  {
311  conn.resize(Quad8::num_nodes);
312  for (auto i : index_range(conn))
313  conn[i] = this->node_id(i);
314 
315  return;
316  }
317 
318  default:
319  libmesh_error_msg("Unsupported IO package " << iop);
320  }
321 }
322 
323 
324 
326 {
327  // This might have curved edges, or might be a curved surface in
328  // 3-space, in which case the full bounding box can be larger than
329  // the bounding box of just the nodes.
330  //
331  //
332  // FIXME - I haven't yet proven the formula below to be correct for
333  // biquadratics - RHS
334  Point pmin, pmax;
335 
336  for (unsigned d=0; d<LIBMESH_DIM; ++d)
337  {
338  Real center = this->point(0)(d);
339  for (unsigned int p=1; p != 8; ++p)
340  center += this->point(p)(d);
341  center /= 8;
342 
343  Real hd = std::abs(center - this->point(0)(d));
344  for (unsigned int p=0; p != 8; ++p)
345  hd = std::max(hd, std::abs(center - this->point(p)(d)));
346 
347  pmin(d) = center - hd;
348  pmax(d) = center + hd;
349  }
350 
351  return BoundingBox(pmin, pmax);
352 }
353 
354 
356 {
357  // This specialization is good for Lagrange mappings only
358  if (this->mapping_type() != LAGRANGE_MAP)
359  return this->Elem::volume();
360 
361  // Make copies of our points. It makes the subsequent calculations a bit
362  // shorter and avoids dereferencing the same pointer multiple times.
363  Point
364  x0 = point(0),
365  x1 = point(1),
366  x2 = point(2),
367  x3 = point(3),
368  x4 = point(4),
369  x5 = point(5),
370  x6 = point(6),
371  x7 = point(7);
372 
373  // Construct constant data vectors.
374  // \vec{x}_{\xi} = \vec{a1}*eta**2 + \vec{b1}*xi*eta + \vec{c1}*xi + \vec{d1}*eta + \vec{e1}
375  // \vec{x}_{\eta} = \vec{a2}*xi**2 + \vec{b2}*xi*eta + \vec{c2}*xi + \vec{d2}*eta + \vec{e2}
376  // This is copy-pasted directly from the output of a Python script.
377  Point
378  a1 = -x0/4 + x1/4 + x2/4 - x3/4 - x5/2 + x7/2,
379  b1 = -x0/2 - x1/2 + x2/2 + x3/2 + x4 - x6,
380  c1 = x0/2 + x1/2 + x2/2 + x3/2 - x4 - x6,
381  d1 = x0/4 - x1/4 + x2/4 - x3/4,
382  e1 = x5/2 - x7/2,
383  a2 = -x0/4 - x1/4 + x2/4 + x3/4 + x4/2 - x6/2,
384  b2 = -x0/2 + x1/2 + x2/2 - x3/2 - x5 + x7,
385  c2 = x0/4 - x1/4 + x2/4 - x3/4,
386  d2 = x0/2 + x1/2 + x2/2 + x3/2 - x5 - x7,
387  e2 = -x4/2 + x6/2;
388 
389  // 3x3 quadrature, exact for bi-quintics
390  const unsigned int N = 3;
391  const Real q[N] = {-std::sqrt(15)/5., 0., std::sqrt(15)/5.};
392  const Real w[N] = {5./9, 8./9, 5./9};
393 
394  Real vol=0.;
395  for (unsigned int i=0; i<N; ++i)
396  for (unsigned int j=0; j<N; ++j)
397  vol += w[i] * w[j] * cross_norm(q[j]*q[j]*a1 + q[i]*q[j]*b1 + q[i]*c1 + q[j]*d1 + e1,
398  q[i]*q[i]*a2 + q[i]*q[j]*b2 + q[i]*c2 + q[j]*d2 + e2);
399 
400  return vol;
401 }
402 
403 
404 
405 unsigned short int Quad8::second_order_adjacent_vertex (const unsigned int n,
406  const unsigned int v) const
407 {
408  libmesh_assert_greater_equal (n, this->n_vertices());
409  libmesh_assert_less (n, this->n_nodes());
410  libmesh_assert_less (v, 2);
411  // use the matrix from \p face_quad.C
412  return _second_order_adjacent_vertices[n-this->n_vertices()][v];
413 }
414 
415 
416 
417 std::pair<unsigned short int, unsigned short int>
418 Quad8::second_order_child_vertex (const unsigned int n) const
419 {
420  libmesh_assert_greater_equal (n, this->n_vertices());
421  libmesh_assert_less (n, this->n_nodes());
422  /*
423  * the _second_order_vertex_child_* vectors are
424  * stored in face_quad.C, since they are identical
425  * for Quad8 and Quad9 (for the first 4 higher-order nodes)
426  */
427  return std::pair<unsigned short int, unsigned short int>
430 }
431 
432 
433 void Quad8::permute(unsigned int perm_num)
434 {
435  libmesh_assert_less (perm_num, 4);
436 
437  for (unsigned int i = 0; i != perm_num; ++i)
438  {
439  swap4nodes(0,1,2,3);
440  swap4nodes(4,5,6,7);
441  swap4neighbors(0,1,2,3);
442  }
443 }
444 
445 
446 void Quad8::flip(BoundaryInfo * boundary_info)
447 {
448  libmesh_assert(boundary_info);
449 
450  swap2nodes(0,1);
451  swap2nodes(2,3);
452  swap2nodes(5,7);
453  swap2neighbors(1,3);
454  swap2boundarysides(1,3,boundary_info);
455  swap2boundaryedges(1,3,boundary_info);
456 }
457 
458 
459 unsigned int Quad8::center_node_on_side(const unsigned short side) const
460 {
461  libmesh_assert_less (side, Quad8::num_sides);
462  return side + 4;
463 }
464 
465 
466 
467 ElemType Quad8::side_type (const unsigned int libmesh_dbg_var(s)) const
468 {
469  libmesh_assert_less (s, 4);
470  return EDGE3;
471 }
472 
473 
474 } // 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 unsigned int n_vertices() const override final
Definition: face_quad.h:102
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
static const int nodes_per_side
Definition: face_quad8.h:191
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Calls cross_norm_sq() and takes the square root of the result.
Definition: type_vector.h:1096
virtual BoundingBox loose_bounding_box() const override
Definition: face_quad8.C:325
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: face_quad8.C:137
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
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
Definition: face_quad8.C:418
virtual dof_id_type key() const override
Definition: face_quad.C:102
virtual unsigned int n_nodes() const override
Definition: face_quad8.h:76
The libMesh namespace provides an interface to certain functionality in the library.
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
Definition: face_quad8.C:446
virtual bool is_vertex(const unsigned int i) const override
Definition: face_quad8.C:108
virtual Order default_order() const override
Definition: face_quad8.C:169
void swap4nodes(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four node_ptrs, "rotating" them.
Definition: elem.h:2143
virtual bool is_face(const unsigned int i) const override
Definition: face_quad8.C:122
ElemMappingType mapping_type() const
Definition: elem.h:3120
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
static const int num_nodes
Geometric constants for Quad8.
Definition: face_quad8.h:190
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: face_quad8.C:238
static const unsigned short int _second_order_vertex_child_index[9]
Vector that names the child vertex index for each second order node.
Definition: face_quad.h:234
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
Definition: face_quad8.C:144
static const unsigned short int _second_order_vertex_child_number[9]
Vector that names a child sharing each second order node.
Definition: face_quad.h:229
static constexpr Real affine_tol
Default tolerance to use in has_affine_map().
Definition: elem.h:2061
static const unsigned short int _second_order_adjacent_vertices[4][2]
Matrix that tells which vertices define the location of mid-side (or second-order) nodes...
Definition: face_quad.h:224
virtual unsigned int n_sub_elem() const override
Definition: face_quad8.h:81
virtual unsigned int n_sides() const override final
Definition: face_quad.h:97
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:2102
virtual bool is_edge(const unsigned int i) const override
Definition: face_quad8.C:115
Defines a Cartesian bounding box by the two corner extremum.
Definition: bounding_box.h:40
static const int num_sides
Geometric constants for Quad4.
Definition: face_quad.h:85
virtual void permute(unsigned int perm_num) override final
Permutes the element (by swapping node and neighbor pointers) according to the specified index...
Definition: face_quad8.C:433
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
Definition: face_quad8.h:197
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: face_quad8.C:127
virtual Real volume() const
Definition: elem.C:3429
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 bool has_affine_map() const override
Definition: face_quad8.C:149
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const override
Definition: face_quad8.C:405
unsigned int center_node_on_side(const unsigned short side) const override final
Definition: face_quad8.C:459
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:3294
virtual Real volume() const override
An optimized method for approximating the area of a QUAD8 using quadrature.
Definition: face_quad8.C:355
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Definition: face_quad8.C:220
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 num_children
Definition: face_quad.h:86
ElemType side_type(const unsigned int s) const override final
Definition: face_quad8.C:467
uint8_t dof_id_type
Definition: id_types.h:67
static const Real _embedding_matrix[num_children][num_nodes][num_nodes]
Matrix that computes new nodal locations/solution values from current nodes/solution.
Definition: face_quad8.h:242
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
Definition: face_quad8.C:209