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