face_quad8.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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/side.h"
20 #include "libmesh/edge_edge3.h"
21 #include "libmesh/face_quad8.h"
23 #include "libmesh/enum_order.h"
24 
25 namespace libMesh
26 {
27 
28 
29 
30 
31 // ------------------------------------------------------------
32 // Quad8 class static member initializations
33 const int Quad8::num_nodes;
34 const int Quad8::num_sides;
35 const int Quad8::num_children;
36 const int Quad8::nodes_per_side;
37 
39  {
40  {0, 1, 4}, // Side 0
41  {1, 2, 5}, // Side 1
42  {2, 3, 6}, // Side 2
43  {3, 0, 7} // Side 3
44  };
45 
46 
47 #ifdef LIBMESH_ENABLE_AMR
48 
50  {
51  // embedding matrix for child 0
52  {
53  // 0 1 2 3 4 5 6 7
54  { 1.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 0
55  { 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000 }, // 1
56  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 2
57  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000 }, // 3
58  { 0.375000, -0.125000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000, 0.00000 }, // 4
59  { -0.187500, -0.187500, -0.187500, -0.187500, 0.750000, 0.375000, 0.250000, 0.375000 }, // 5
60  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.250000, 0.375000, 0.750000 }, // 6
61  { 0.375000, 0.00000, 0.00000, -0.125000, 0.00000, 0.00000, 0.00000, 0.750000 } // 7
62  },
63 
64  // embedding matrix for child 1
65  {
66  // 0 1 2 3 4 5 6 7
67  { 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000 }, // 0
68  { 0.00000, 1.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 1
69  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000 }, // 2
70  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 3
71  { -0.125000, 0.375000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000, 0.00000 }, // 4
72  { 0.00000, 0.375000, -0.125000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000 }, // 5
73  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.750000, 0.375000, 0.250000 }, // 6
74  { -0.187500, -0.187500, -0.187500, -0.187500, 0.750000, 0.375000, 0.250000, 0.375000 } // 7
75  },
76 
77  // embedding matrix for child 2
78  {
79  // 0 1 2 3 4 5 6 7
80  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000 }, // 0
81  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 1
82  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000 }, // 2
83  { 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 3
84  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.250000, 0.375000, 0.750000 }, // 4
85  { -0.187500, -0.187500, -0.187500, -0.187500, 0.250000, 0.375000, 0.750000, 0.375000 }, // 5
86  { 0.00000, 0.00000, -0.125000, 0.375000, 0.00000, 0.00000, 0.750000, 0.00000 }, // 6
87  { -0.125000, 0.00000, 0.00000, 0.375000, 0.00000, 0.00000, 0.00000, 0.750000 } // 7
88  },
89 
90  // embedding matrix for child 3
91  {
92  // 0 1 2 3 4 5 6 7
93  { -0.250000, -0.250000, -0.250000, -0.250000, 0.500000, 0.500000, 0.500000, 0.500000 }, // 0
94  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000, 0.00000 }, // 1
95  { 0.00000, 0.00000, 1.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000 }, // 2
96  { 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 0.00000, 1.00000, 0.00000 }, // 3
97  { -0.187500, -0.187500, -0.187500, -0.187500, 0.375000, 0.750000, 0.375000, 0.250000 }, // 4
98  { 0.00000, -0.125000, 0.375000, 0.00000, 0.00000, 0.750000, 0.00000, 0.00000 }, // 5
99  { 0.00000, 0.00000, 0.375000, -0.125000, 0.00000, 0.00000, 0.750000, 0.00000 }, // 6
100  { -0.187500, -0.187500, -0.187500, -0.187500, 0.250000, 0.375000, 0.750000, 0.375000 } // 7
101  }
102  };
103 
104 
105 #endif
106 
107 
108 // ------------------------------------------------------------
109 // Quad8 class member functions
110 
111 bool Quad8::is_vertex(const unsigned int i) const
112 {
113  if (i < 4)
114  return true;
115  return false;
116 }
117 
118 bool Quad8::is_edge(const unsigned int i) const
119 {
120  if (i < 4)
121  return false;
122  return true;
123 }
124 
125 bool Quad8::is_face(const unsigned int) const
126 {
127  return false;
128 }
129 
130 bool Quad8::is_node_on_side(const unsigned int n,
131  const unsigned int s) const
132 {
133  libmesh_assert_less (s, n_sides());
134  return std::find(std::begin(side_nodes_map[s]),
136  n) != std::end(side_nodes_map[s]);
137 }
138 
139 std::vector<unsigned>
140 Quad8::nodes_on_side(const unsigned int s) const
141 {
142  libmesh_assert_less(s, n_sides());
143  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s])};
144 }
145 
147 {
148  // make sure corners form a parallelogram
149  Point v = this->point(1) - this->point(0);
150  if (!v.relative_fuzzy_equals(this->point(2) - this->point(3)))
151  return false;
152  // make sure sides are straight
153  v /= 2;
154  if (!v.relative_fuzzy_equals(this->point(4) - this->point(0)) ||
155  !v.relative_fuzzy_equals(this->point(6) - this->point(3)))
156  return false;
157  v = (this->point(3) - this->point(0))/2;
158  if (!v.relative_fuzzy_equals(this->point(7) - this->point(0)) ||
159  !v.relative_fuzzy_equals(this->point(5) - this->point(1)))
160  return false;
161  return true;
162 }
163 
164 
165 
167 {
168  return SECOND;
169 }
170 
171 
172 
173 dof_id_type Quad8::key (const unsigned int s) const
174 {
175  libmesh_assert_less (s, this->n_sides());
176 
177  switch (s)
178  {
179  case 0:
180 
181  return
182  this->compute_key (this->node_id(4));
183 
184  case 1:
185 
186  return
187  this->compute_key (this->node_id(5));
188 
189  case 2:
190 
191  return
192  this->compute_key (this->node_id(6));
193 
194  case 3:
195 
196  return
197  this->compute_key (this->node_id(7));
198 
199  default:
200  libmesh_error_msg("Invalid side s = " << s);
201  }
202 }
203 
204 
205 
206 unsigned int Quad8::which_node_am_i(unsigned int side,
207  unsigned int side_node) const
208 {
209  libmesh_assert_less (side, this->n_sides());
210  libmesh_assert_less (side_node, Quad8::nodes_per_side);
211 
212  return Quad8::side_nodes_map[side][side_node];
213 }
214 
215 
216 
217 std::unique_ptr<Elem> Quad8::build_side_ptr (const unsigned int i,
218  bool proxy)
219 {
220  libmesh_assert_less (i, this->n_sides());
221 
222  if (proxy)
223  return libmesh_make_unique<Side<Edge3,Quad8>>(this,i);
224 
225  else
226  {
227  std::unique_ptr<Elem> edge = libmesh_make_unique<Edge3>();
228  edge->subdomain_id() = this->subdomain_id();
229 
230  // Set the nodes
231  for (unsigned n=0; n<edge->n_nodes(); ++n)
232  edge->set_node(n) = this->node_ptr(Quad8::side_nodes_map[i][n]);
233 
234  return edge;
235  }
236 }
237 
238 
239 
240 void Quad8::build_side_ptr (std::unique_ptr<Elem> & side,
241  const unsigned int i)
242 {
243  this->simple_build_side_ptr<Quad8>(side, i, EDGE3);
244 }
245 
246 
247 
248 
249 
250 
251 void Quad8::connectivity(const unsigned int sf,
252  const IOPackage iop,
253  std::vector<dof_id_type> & conn) const
254 {
255  libmesh_assert_less (sf, this->n_sub_elem());
256  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
257 
258  switch (iop)
259  {
260  // Note: TECPLOT connectivity is output as four triangles with
261  // a central quadrilateral. Therefore, the first four connectivity
262  // arrays are degenerate quads (triangles in Tecplot).
263  case TECPLOT:
264  {
265  // Create storage
266  conn.resize(4);
267 
268  switch(sf)
269  {
270  case 0:
271  // linear sub-tri 0
272  conn[0] = this->node_id(0)+1;
273  conn[1] = this->node_id(4)+1;
274  conn[2] = this->node_id(7)+1;
275  conn[3] = this->node_id(7)+1;
276 
277  return;
278 
279  case 1:
280  // linear sub-tri 1
281  conn[0] = this->node_id(4)+1;
282  conn[1] = this->node_id(1)+1;
283  conn[2] = this->node_id(5)+1;
284  conn[3] = this->node_id(5)+1;
285 
286  return;
287 
288  case 2:
289  // linear sub-tri 2
290  conn[0] = this->node_id(5)+1;
291  conn[1] = this->node_id(2)+1;
292  conn[2] = this->node_id(6)+1;
293  conn[3] = this->node_id(6)+1;
294 
295  return;
296 
297  case 3:
298  // linear sub-tri 3
299  conn[0] = this->node_id(7)+1;
300  conn[1] = this->node_id(6)+1;
301  conn[2] = this->node_id(3)+1;
302  conn[3] = this->node_id(3)+1;
303 
304  return;
305 
306  case 4:
307  // linear sub-quad
308  conn[0] = this->node_id(4)+1;
309  conn[1] = this->node_id(5)+1;
310  conn[2] = this->node_id(6)+1;
311  conn[3] = this->node_id(7)+1;
312 
313  return;
314 
315  default:
316  libmesh_error_msg("Invalid sf = " << sf);
317  }
318  }
319 
320 
321  // Note: VTK connectivity is output as four triangles with
322  // a central quadrilateral. Therefore most of the connectivity
323  // arrays have length three.
324  case VTK:
325  {
326  // Create storage
327  conn.resize(8);
328  conn[0] = this->node_id(0);
329  conn[1] = this->node_id(1);
330  conn[2] = this->node_id(2);
331  conn[3] = this->node_id(3);
332  conn[4] = this->node_id(4);
333  conn[5] = this->node_id(5);
334  conn[6] = this->node_id(6);
335  conn[7] = this->node_id(7);
336  return;
337  /*
338  conn.resize(3);
339 
340  switch (sf)
341  {
342  case 0:
343  // linear sub-tri 0
344  conn[0] = this->node_id(0);
345  conn[1] = this->node_id(4);
346  conn[2] = this->node_id(7);
347 
348  return;
349 
350  case 1:
351  // linear sub-tri 1
352  conn[0] = this->node_id(4);
353  conn[1] = this->node_id(1);
354  conn[2] = this->node_id(5);
355 
356  return;
357 
358  case 2:
359  // linear sub-tri 2
360  conn[0] = this->node_id(5);
361  conn[1] = this->node_id(2);
362  conn[2] = this->node_id(6);
363 
364  return;
365 
366  case 3:
367  // linear sub-tri 3
368  conn[0] = this->node_id(7);
369  conn[1] = this->node_id(6);
370  conn[2] = this->node_id(3);
371 
372  return;
373 
374  case 4:
375  conn.resize(4);
376 
377  // linear sub-quad
378  conn[0] = this->node_id(4);
379  conn[1] = this->node_id(5);
380  conn[2] = this->node_id(6);
381  conn[3] = this->node_id(7);
382  */
383  // return;
384 
385  // default:
386  // libmesh_error_msg("Invalid sf = " << sf);
387  // }
388  }
389 
390  default:
391  libmesh_error_msg("Unsupported IO package " << iop);
392  }
393 }
394 
395 
396 
398 {
399  // This might have curved edges, or might be a curved surface in
400  // 3-space, in which case the full bounding box can be larger than
401  // the bounding box of just the nodes.
402  //
403  //
404  // FIXME - I haven't yet proven the formula below to be correct for
405  // biquadratics - RHS
406  Point pmin, pmax;
407 
408  for (unsigned d=0; d<LIBMESH_DIM; ++d)
409  {
410  Real center = this->point(0)(d);
411  for (unsigned int p=1; p != 8; ++p)
412  center += this->point(p)(d);
413  center /= 8;
414 
415  Real hd = std::abs(center - this->point(0)(d));
416  for (unsigned int p=0; p != 8; ++p)
417  hd = std::max(hd, std::abs(center - this->point(p)(d)));
418 
419  pmin(d) = center - hd;
420  pmax(d) = center + hd;
421  }
422 
423  return BoundingBox(pmin, pmax);
424 }
425 
426 
428 {
429  // Make copies of our points. It makes the subsequent calculations a bit
430  // shorter and avoids dereferencing the same pointer multiple times.
431  Point
432  x0 = point(0),
433  x1 = point(1),
434  x2 = point(2),
435  x3 = point(3),
436  x4 = point(4),
437  x5 = point(5),
438  x6 = point(6),
439  x7 = point(7);
440 
441  // Construct constant data vectors.
442  // \vec{x}_{\xi} = \vec{a1}*eta**2 + \vec{b1}*xi*eta + \vec{c1}*xi + \vec{d1}*eta + \vec{e1}
443  // \vec{x}_{\eta} = \vec{a2}*xi**2 + \vec{b2}*xi*eta + \vec{c2}*xi + \vec{d2}*eta + \vec{e2}
444  // This is copy-pasted directly from the output of a Python script.
445  Point
446  a1 = -x0/4 + x1/4 + x2/4 - x3/4 - x5/2 + x7/2,
447  b1 = -x0/2 - x1/2 + x2/2 + x3/2 + x4 - x6,
448  c1 = x0/2 + x1/2 + x2/2 + x3/2 - x4 - x6,
449  d1 = x0/4 - x1/4 + x2/4 - x3/4,
450  e1 = x5/2 - x7/2,
451  a2 = -x0/4 - x1/4 + x2/4 + x3/4 + x4/2 - x6/2,
452  b2 = -x0/2 + x1/2 + x2/2 - x3/2 - x5 + x7,
453  c2 = x0/4 - x1/4 + x2/4 - x3/4,
454  d2 = x0/2 + x1/2 + x2/2 + x3/2 - x5 - x7,
455  e2 = -x4/2 + x6/2;
456 
457  // 3x3 quadrature, exact for bi-quintics
458  const unsigned int N = 3;
459  const Real q[N] = {-std::sqrt(15)/5., 0., std::sqrt(15)/5.};
460  const Real w[N] = {5./9, 8./9, 5./9};
461 
462  Real vol=0.;
463  for (unsigned int i=0; i<N; ++i)
464  for (unsigned int j=0; j<N; ++j)
465  vol += w[i] * w[j] * cross_norm(q[j]*q[j]*a1 + q[i]*q[j]*b1 + q[i]*c1 + q[j]*d1 + e1,
466  q[i]*q[i]*a2 + q[i]*q[j]*b2 + q[i]*c2 + q[j]*d2 + e2);
467 
468  return vol;
469 }
470 
471 
472 
473 unsigned short int Quad8::second_order_adjacent_vertex (const unsigned int n,
474  const unsigned int v) const
475 {
476  libmesh_assert_greater_equal (n, this->n_vertices());
477  libmesh_assert_less (n, this->n_nodes());
478  libmesh_assert_less (v, 2);
479  // use the matrix from \p face_quad.C
480  return _second_order_adjacent_vertices[n-this->n_vertices()][v];
481 }
482 
483 
484 
485 std::pair<unsigned short int, unsigned short int>
486 Quad8::second_order_child_vertex (const unsigned int n) const
487 {
488  libmesh_assert_greater_equal (n, this->n_vertices());
489  libmesh_assert_less (n, this->n_nodes());
490  /*
491  * the _second_order_vertex_child_* vectors are
492  * stored in face_quad.C, since they are identical
493  * for Quad8 and Quad9 (for the first 4 higher-order nodes)
494  */
495  return std::pair<unsigned short int, unsigned short int>
498 }
499 
500 } // namespace libMesh
double abs(double a)
virtual unsigned int n_vertices() const override final
Definition: face_quad.h:96
static const int num_sides
Definition: face_quad8.h:188
static const int nodes_per_side
Definition: face_quad8.h:190
T cross_norm(const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1097
virtual unsigned int which_node_am_i(unsigned int side, unsigned int side_node) const override
Definition: face_quad8.C:206
virtual BoundingBox loose_bounding_box() const override
Definition: face_quad8.C:397
unsigned short int side
Definition: xdr_io.C:50
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
Definition: face_quad8.C:140
virtual std::pair< unsigned short int, unsigned short int > second_order_child_vertex(const unsigned int n) const override
Definition: face_quad8.C:486
virtual dof_id_type key() const override
Definition: face_quad.C:76
IterBase * end
virtual unsigned int n_nodes() const override
Definition: face_quad8.h:76
long double max(long double a, double b)
virtual bool is_vertex(const unsigned int i) const override
Definition: face_quad8.C:111
virtual Order default_order() const override
Definition: face_quad8.C:166
static const float _embedding_matrix[num_children][num_nodes][num_nodes]
Definition: face_quad8.h:233
virtual bool is_face(const unsigned int i) const override
Definition: face_quad8.C:125
static const int num_nodes
Definition: face_quad8.h:187
virtual void connectivity(const unsigned int sf, const IOPackage iop, std::vector< dof_id_type > &conn) const override
Definition: face_quad8.C:251
static const unsigned short int _second_order_vertex_child_index[9]
Definition: face_quad.h:202
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy) override
Definition: face_quad8.C:217
static const unsigned short int _second_order_vertex_child_number[9]
Definition: face_quad.h:197
static const unsigned short int _second_order_adjacent_vertices[4][2]
Definition: face_quad.h:192
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:91
virtual bool is_edge(const unsigned int i) const override
Definition: face_quad8.C:118
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
subdomain_id_type subdomain_id() const
Definition: elem.h:2034
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
Definition: face_quad8.h:196
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
Definition: face_quad8.C:130
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:1957
virtual bool has_affine_map() const override
Definition: face_quad8.C:146
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const override
Definition: face_quad8.C:473
static const int num_children
Definition: face_quad8.h:189
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:2754
virtual Real volume() const override
Definition: face_quad8.C:427
A geometric point in (x,y,z) space.
Definition: point.h:38
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:1914
const Point & point(const unsigned int i) const
Definition: elem.h:1892
bool relative_fuzzy_equals(const TypeVector< T > &rhs, Real tol=TOLERANCE) const
Definition: type_vector.h:990
std::unique_ptr< Elem > side(const unsigned int i) const
Definition: elem.h:2202
uint8_t dof_id_type
Definition: id_types.h:64