face_quad.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/face_quad.h"
20 #include "libmesh/edge_edge2.h"
21 #include "libmesh/face_quad4.h"
23 
24 // C++ includes
25 #include <array>
26 
27 
28 namespace libMesh
29 {
30 
31 
32 
33 // ------------------------------------------------------------
34 // Quad class static member initializations
35 
36 
37 // We need to require C++11...
38 const Real Quad::_master_points[9][3] =
39  {
40  {-1, -1},
41  {1, -1},
42  {1, 1},
43  {-1, 1},
44  {0, -1},
45  {1, 0},
46  {0, 1},
47  {-1, 0},
48  {0, 0}
49  };
50 
51 
52 
53 // ------------------------------------------------------------
54 // Quad class member functions
55 dof_id_type Quad::key (const unsigned int s) const
56 {
57  libmesh_assert_less (s, this->n_sides());
58 
59  return this->compute_key(this->node_id(Quad4::side_nodes_map[s][0]),
60  this->node_id(Quad4::side_nodes_map[s][1]));
61 }
62 
63 
64 
65 unsigned int Quad::which_node_am_i(unsigned int side,
66  unsigned int side_node) const
67 {
68  libmesh_assert_less (side, this->n_sides());
69  libmesh_assert_less (side_node, 2);
70 
71  return Quad4::side_nodes_map[side][side_node];
72 }
73 
74 
75 
77 {
78  return this->compute_key(this->node_id(0),
79  this->node_id(1),
80  this->node_id(2),
81  this->node_id(3));
82 }
83 
84 
85 
86 std::unique_ptr<Elem> Quad::side_ptr (const unsigned int i)
87 {
88  libmesh_assert_less (i, this->n_sides());
89 
90  std::unique_ptr<Elem> edge = libmesh_make_unique<Edge2>();
91 
92  for (unsigned n=0; n<edge->n_nodes(); ++n)
93  edge->set_node(n) = this->node_ptr(Quad4::side_nodes_map[i][n]);
94 
95  return edge;
96 }
97 
98 
99 
100 void Quad::side_ptr (std::unique_ptr<Elem> & side,
101  const unsigned int i)
102 {
103  this->simple_build_side_ptr<Quad4>(side, i, EDGE2);
104 }
105 
106 
107 
108 bool Quad::is_child_on_side(const unsigned int c,
109  const unsigned int s) const
110 {
111  libmesh_assert_less (c, this->n_children());
112  libmesh_assert_less (s, this->n_sides());
113 
114  // A quad's children and nodes don't share the same ordering:
115  // child 2 and 3 are swapped;
116  unsigned int n = (c < 2) ? c : 5-c;
117  return (n == s || n == (s+1)%4);
118 }
119 
120 
121 
122 unsigned int Quad::opposite_side(const unsigned int side_in) const
123 {
124  libmesh_assert_less (side_in, 4);
125 
126  return (side_in + 2) % 4;
127 }
128 
129 
130 
131 unsigned int Quad::opposite_node(const unsigned int node_in,
132  const unsigned int side_in) const
133 {
134  libmesh_assert_less (node_in, 8);
135  libmesh_assert_less (node_in, this->n_nodes());
136  libmesh_assert_less (side_in, this->n_sides());
137  libmesh_assert(this->is_node_on_side(node_in, side_in));
138 
139  static const unsigned char side02_nodes_map[] =
140  {3, 2, 1, 0, 6, 255, 4, 255};
141  static const unsigned char side13_nodes_map[] =
142  {1, 0, 3, 2, 255, 7, 255, 5};
143 
144  switch (side_in)
145  {
146  case 0:
147  case 2:
148  return side02_nodes_map[node_in];
149  case 1:
150  case 3:
151  return side13_nodes_map[node_in];
152  default:
153  libmesh_error_msg("Unsupported side_in = " << side_in);
154  }
155 }
156 
157 
159 {
160  switch (q)
161  {
162  // CUBIT 15.1 User Documentation:
163  // Aspect Ratio: Maximum edge length ratios
164  case ASPECT_RATIO:
165  {
166  Real lengths[4] = {this->length(0,1), this->length(1,2), this->length(2,3), this->length(3,0)};
167  Real
168  max = *std::max_element(lengths, lengths+4),
169  min = *std::min_element(lengths, lengths+4);
170 
171  // Return 0. instead of dividing by zero.
172  if (min == 0.)
173  return 0.;
174  else
175  return max / min;
176  }
177 
178  // Compute the min/max diagonal ratio.
179  // This is modeled after the Hex element
180  case DISTORTION:
181  case DIAGONAL:
182  {
183  // Diagonal between node 0 and node 2
184  const Real d02 = this->length(0,2);
185 
186  // Diagonal between node 1 and node 3
187  const Real d13 = this->length(1,3);
188 
189  // Find the biggest and smallest diagonals
190  if ((d02 > 0.) && (d13 >0.))
191  if (d02 < d13) return d02 / d13;
192  else return d13 / d02;
193  else
194  return 0.;
195  break;
196  }
197 
198  // CUBIT 15.1 User Documentation:
199  // Stretch: Sqrt(2) * minimum edge length / maximum diagonal length
200  case STRETCH:
201  {
202  Real lengths[4] = {this->length(0,1), this->length(1,2), this->length(2,3), this->length(3,0)};
203  Real min_edge = *std::min_element(lengths, lengths+4);
204  Real d_max = std::max(this->length(0,2), this->length(1,3));
205 
206  // Return 0. instead of dividing by zero.
207  if (d_max == 0.)
208  return 0.;
209  else
210  return std::sqrt(2) * min_edge / d_max;
211  }
212 
213  case SHAPE:
214  case SKEW:
215  {
216  // From: P. Knupp, "Algebraic mesh quality metrics for
217  // unstructured initial meshes," Finite Elements in Analysis
218  // and Design 39, 2003, p. 217-241, Sections 5.2 and 5.3.
219  typedef std::array<Real, 4> Array4;
220  typedef std::array<Real, 6> Array6;
221 
222  // x, y, z node coordinates.
223  std::vector<Real>
224  x = {this->point(0)(0), this->point(1)(0), this->point(2)(0), this->point(3)(0)},
225  y = {this->point(0)(1), this->point(1)(1), this->point(2)(1), this->point(3)(1)},
226  z = {this->point(0)(2), this->point(1)(2), this->point(2)(2), this->point(3)(2)};
227 
228  // Nodal Jacobians. These are 3x2 matrices, hence we represent
229  // them by Array6.
230  std::array<Array6, 4> A;
231  for (unsigned int k=0; k<4; ++k)
232  {
233  unsigned int
234  kp1 = k+1 > 3 ? k+1-4 : k+1,
235  kp3 = k+3 > 3 ? k+3-4 : k+3;
236 
237  // To initialize std::array we need double curly braces in
238  // C++11 but not C++14 apparently.
239  A[k] = {{x[kp1] - x[k], x[kp3] - x[k],
240  y[kp1] - y[k], y[kp3] - y[k],
241  z[kp1] - z[k], z[kp3] - z[k]}};
242  }
243 
244  // Compute metric tensors, T_k = A_k^T * A_k. These are 2x2
245  // square matrices, hence we represent them by Array4.
246  std::array<Array4, 4> T;
247  for (unsigned int k=0; k<4; ++k)
248  {
249  Real
250  top_left = A[k][0]*A[k][0] + A[k][2]*A[k][2] + A[k][4]*A[k][4],
251  off_diag = A[k][0]*A[k][1] + A[k][2]*A[k][3] + A[k][4]*A[k][5],
252  bot_rigt = A[k][1]*A[k][1] + A[k][3]*A[k][3] + A[k][5]*A[k][5];
253 
254  T[k] = {{top_left, off_diag,
255  off_diag, bot_rigt}};
256  }
257 
258 
259  // Nodal areas. These are approximated as sqrt(det(A^T * A))
260  // to handle the general case of a 2D element living in 3D.
261  Array4 alpha;
262  for (unsigned int k=0; k<4; ++k)
263  alpha[k] = std::sqrt(T[k][0]*T[k][3] - T[k][1]*T[k][2]);
264 
265  // All nodal areas must be strictly positive. Return 0 (the
266  // lowest quality) otherwise. We know they can't be negative
267  // because they are the result of a sqrt, but they might be
268  // identically 0.
269  if (*std::min_element(alpha.begin(), alpha.end()) == 0.)
270  return 0.;
271 
272  // Compute and return the shape metric. These only use the
273  // diagonal entries of the T_k.
274  Real den = 0.;
275  if (q == SHAPE)
276  {
277  for (unsigned int k=0; k<4; ++k)
278  den += (T[k][0] + T[k][3]) / alpha[k];
279  return (den == 0.) ? 0 : (8. / den);
280  }
281  else
282  {
283  for (unsigned int k=0; k<4; ++k)
284  den += std::sqrt(T[k][0] * T[k][3]) / alpha[k];
285  return (den == 0.) ? 0 : (4. / den);
286  }
287  }
288 
289  default:
290  return Elem::quality(q);
291  }
292 }
293 
294 
295 
296 
297 
298 
299 std::pair<Real, Real> Quad::qual_bounds (const ElemQuality q) const
300 {
301  std::pair<Real, Real> bounds;
302 
303  switch (q)
304  {
305 
306  case ASPECT_RATIO:
307  bounds.first = 1.;
308  bounds.second = 4.;
309  break;
310 
311  case SKEW:
312  bounds.first = 0.;
313  bounds.second = 0.5;
314  break;
315 
316  case TAPER:
317  bounds.first = 0.;
318  bounds.second = 0.7;
319  break;
320 
321  case WARP:
322  bounds.first = 0.9;
323  bounds.second = 1.;
324  break;
325 
326  case STRETCH:
327  bounds.first = 0.25;
328  bounds.second = 1.;
329  break;
330 
331  case MIN_ANGLE:
332  bounds.first = 45.;
333  bounds.second = 90.;
334  break;
335 
336  case MAX_ANGLE:
337  bounds.first = 90.;
338  bounds.second = 135.;
339  break;
340 
341  case CONDITION:
342  bounds.first = 1.;
343  bounds.second = 4.;
344  break;
345 
346  case JACOBIAN:
347  bounds.first = 0.5;
348  bounds.second = 1.;
349  break;
350 
351  case SHEAR:
352  case SHAPE:
353  case SIZE:
354  bounds.first = 0.3;
355  bounds.second = 1.;
356  break;
357 
358  case DISTORTION:
359  bounds.first = 0.6;
360  bounds.second = 1.;
361  break;
362 
363  default:
364  libMesh::out << "Warning: Invalid quality measure chosen." << std::endl;
365  bounds.first = -1;
366  bounds.second = -1;
367  }
368 
369  return bounds;
370 }
371 
372 
373 
374 
375 const unsigned short int Quad::_second_order_adjacent_vertices[4][2] =
376  {
377  {0, 1}, // vertices adjacent to node 4
378  {1, 2}, // vertices adjacent to node 5
379  {2, 3}, // vertices adjacent to node 6
380  {0, 3} // vertices adjacent to node 7
381  };
382 
383 
384 
385 const unsigned short int Quad::_second_order_vertex_child_number[9] =
386  {
387  99,99,99,99, // Vertices
388  0,1,2,0, // Edges
389  0 // Interior
390  };
391 
392 
393 
394 const unsigned short int Quad::_second_order_vertex_child_index[9] =
395  {
396  99,99,99,99, // Vertices
397  1,2,3,3, // Edges
398  2 // Interior
399  };
400 
401 
402 
403 #ifdef LIBMESH_ENABLE_AMR
404 
405 // We number 25 "possible node locations" for a 2x2 refinement of
406 // quads with up to 3x3 nodes each
407 const int Quad::_child_node_lookup[4][9] =
408  {
409  // node lookup for child 0 (near node 0)
410  { 0, 2, 12, 10, 1, 7, 11, 5, 6},
411 
412  // node lookup for child 1 (near node 1)
413  { 2, 4, 14, 12, 3, 9, 13, 7, 8},
414 
415  // node lookup for child 2 (near node 3)
416  { 10, 12, 22, 20, 11, 17, 21, 15, 16},
417 
418  // node lookup for child 3 (near node 2)
419  { 12, 14, 24, 22, 13, 19, 23, 17, 18}
420  };
421 
422 
423 #endif
424 
425 } // namespace libMesh
virtual std::pair< Real, Real > qual_bounds(const ElemQuality q) const override
Definition: face_quad.C:299
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
Definition: face_quad4.h:146
static const int _child_node_lookup[4][9]
Definition: face_quad.h:213
unsigned short int side
Definition: xdr_io.C:50
virtual std::unique_ptr< Elem > side_ptr(const unsigned int i) override final
Definition: face_quad.C:86
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const =0
virtual dof_id_type key() const override
Definition: face_quad.C:76
static const Real _master_points[9][3]
Definition: face_quad.h:207
long double max(long double a, double b)
virtual unsigned int n_children() const override final
Definition: face_quad.h:106
Real length(const unsigned int n1, const unsigned int n2) const
Definition: elem.C:390
static const unsigned short int _second_order_vertex_child_index[9]
Definition: face_quad.h:202
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 Real quality(const ElemQuality q) const override
Definition: face_quad.C:158
virtual unsigned int n_sides() const override final
Definition: face_quad.h:91
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned int which_node_am_i(unsigned int side, unsigned int side_node) const override
Definition: face_quad.C:65
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:1957
static PetscErrorCode Mat * A
virtual Real quality(const ElemQuality q) const
Definition: elem.C:1403
virtual unsigned int opposite_side(const unsigned int s) const override final
Definition: face_quad.C:122
virtual unsigned int opposite_node(const unsigned int n, const unsigned int s) const override final
Definition: face_quad.C:131
virtual unsigned int n_nodes() const override
Definition: face_quad.h:86
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:2754
OStreamProxy out(std::cout)
long double min(long double a, double b)
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
std::unique_ptr< Elem > side(const unsigned int i) const
Definition: elem.h:2202
virtual bool is_child_on_side(const unsigned int c, const unsigned int s) const override final
Definition: face_quad.C:108
uint8_t dof_id_type
Definition: id_types.h:64