libMesh
cell_pyramid14.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 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 // C++ includes
20 
21 // Local includes
22 #include "libmesh/side.h"
23 #include "libmesh/cell_pyramid14.h"
24 #include "libmesh/edge_edge3.h"
25 #include "libmesh/face_tri6.h"
26 #include "libmesh/face_quad9.h"
27 
28 namespace libMesh
29 {
30 
31 
32 
33 
34 // ------------------------------------------------------------
35 // Pyramid14 class static member initializations
36 const unsigned int Pyramid14::side_nodes_map[5][9] =
37  {
38  {0, 1, 4, 5, 10, 9, 99, 99, 99}, // Side 0 (front)
39  {1, 2, 4, 6, 11, 10, 99, 99, 99}, // Side 1 (right)
40  {2, 3, 4, 7, 12, 11, 99, 99, 99}, // Side 2 (back)
41  {3, 0, 4, 8, 9, 12, 99, 99, 99}, // Side 3 (left)
42  {0, 3, 2, 1, 8, 7, 6, 5, 13} // Side 4 (base)
43  };
44 
45 const unsigned int Pyramid14::edge_nodes_map[8][3] =
46  {
47  {0, 1, 5}, // Edge 0
48  {1, 2, 6}, // Edge 1
49  {2, 3, 7}, // Edge 2
50  {0, 3, 8}, // Edge 3
51  {0, 4, 9}, // Edge 4
52  {1, 4, 10}, // Edge 5
53  {2, 4, 11}, // Edge 6
54  {3, 4, 12} // Edge 7
55  };
56 
57 
58 
59 // ------------------------------------------------------------
60 // Pyramid14 class member functions
61 
62 bool Pyramid14::is_vertex(const unsigned int i) const
63 {
64  if (i < 5)
65  return true;
66  return false;
67 }
68 
69 
70 
71 bool Pyramid14::is_edge(const unsigned int i) const
72 {
73  if (i < 5)
74  return false;
75  if (i == 13)
76  return false;
77  return true;
78 }
79 
80 
81 
82 bool Pyramid14::is_face(const unsigned int i) const
83 {
84  if (i == 13)
85  return true;
86  return false;
87 }
88 
89 
90 
91 bool Pyramid14::is_node_on_side(const unsigned int n,
92  const unsigned int s) const
93 {
94  libmesh_assert_less (s, n_sides());
95  for (unsigned int i = 0; i != 9; ++i)
96  if (side_nodes_map[s][i] == n)
97  return true;
98  return false;
99 }
100 
101 bool Pyramid14::is_node_on_edge(const unsigned int n,
102  const unsigned int e) const
103 {
104  libmesh_assert_less (e, n_edges());
105  for (unsigned int i = 0; i != 3; ++i)
106  if (edge_nodes_map[e][i] == n)
107  return true;
108  return false;
109 }
110 
111 
112 
113 bool Pyramid14::has_affine_map() const
114 {
115  // TODO: If the base is a parallelogram and all the triangular faces are planar,
116  // the map should be linear, but I need to test this theory...
117  return false;
118 }
119 
120 
121 
122 dof_id_type Pyramid14::key (const unsigned int s) const
123 {
124  libmesh_assert_less (s, this->n_sides());
125 
126  switch (s)
127  {
128  case 0: // triangular face 1
129  case 1: // triangular face 2
130  case 2: // triangular face 3
131  case 3: // triangular face 4
132  return Pyramid::key(s);
133 
134  case 4: // the quad face at z=0
135  return this->compute_key (this->node_id(13));
136 
137  default:
138  libmesh_error_msg("Invalid side s = " << s);
139  }
140 
141  libmesh_error_msg("We'll never get here!");
142  return 0;
143 }
144 
145 
146 
147 unsigned int Pyramid14::which_node_am_i(unsigned int side,
148  unsigned int side_node) const
149 {
150  libmesh_assert_less (side, this->n_sides());
151 
152  // Never more than 9 nodes per side.
153  libmesh_assert_less(side_node, 9);
154 
155  // Some sides have 6 nodes.
156  libmesh_assert(side == 4 || side_node < 6);
157 
158  return Pyramid14::side_nodes_map[side][side_node];
159 }
160 
161 
162 
163 UniquePtr<Elem> Pyramid14::build_side_ptr (const unsigned int i, bool proxy)
164 {
165  libmesh_assert_less (i, this->n_sides());
166 
167  if (proxy)
168  {
169  switch (i)
170  {
171  case 0:
172  case 1:
173  case 2:
174  case 3:
175  return UniquePtr<Elem>(new Side<Tri6,Pyramid14>(this,i));
176 
177  case 4:
178  return UniquePtr<Elem>(new Side<Quad9,Pyramid14>(this,i));
179 
180  default:
181  libmesh_error_msg("Invalid side i = " << i);
182  }
183  }
184 
185  else
186  {
187  // Create NULL pointer to be initialized, returned later.
188  Elem * face = libmesh_nullptr;
189 
190  switch (i)
191  {
192  case 0: // triangular face 1
193  case 1: // triangular face 2
194  case 2: // triangular face 3
195  case 3: // triangular face 4
196  {
197  face = new Tri6;
198  break;
199  }
200  case 4: // the quad face at z=0
201  {
202  face = new Quad9;
203  break;
204  }
205  default:
206  libmesh_error_msg("Invalid side i = " << i);
207  }
208 
209  face->subdomain_id() = this->subdomain_id();
210 
211  // Set the nodes
212  for (unsigned n=0; n<face->n_nodes(); ++n)
213  face->set_node(n) = this->node_ptr(Pyramid14::side_nodes_map[i][n]);
214 
215  return UniquePtr<Elem>(face);
216  }
217 
218  libmesh_error_msg("We'll never get here!");
219  return UniquePtr<Elem>();
220 }
221 
222 
223 
224 UniquePtr<Elem> Pyramid14::build_edge_ptr (const unsigned int i)
225 {
226  libmesh_assert_less (i, this->n_edges());
227 
228  return UniquePtr<Elem>(new SideEdge<Edge3,Pyramid14>(this,i));
229 }
230 
231 
232 
233 void Pyramid14::connectivity(const unsigned int libmesh_dbg_var(sc),
234  const IOPackage iop,
235  std::vector<dof_id_type> & /*conn*/) const
236 {
237  libmesh_assert(_nodes);
238  libmesh_assert_less (sc, this->n_sub_elem());
239  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
240 
241  switch (iop)
242  {
243  case TECPLOT:
244  {
245  // TODO
246  libmesh_not_implemented();
247  }
248 
249  case VTK:
250  {
251  // TODO
252  libmesh_not_implemented();
253  }
254 
255  default:
256  libmesh_error_msg("Unsupported IO package " << iop);
257  }
258 }
259 
260 
261 
262 unsigned int Pyramid14::n_second_order_adjacent_vertices (const unsigned int n) const
263 {
264  switch (n)
265  {
266  case 5:
267  case 6:
268  case 7:
269  case 8:
270  case 9:
271  case 10:
272  case 11:
273  case 12:
274  return 2;
275 
276  case 13:
277  return 4;
278 
279  default:
280  libmesh_error_msg("Invalid node n = " << n);
281  }
282 
283  libmesh_error_msg("We'll never get here!");
284  return libMesh::invalid_uint;
285 }
286 
287 
288 unsigned short int Pyramid14::second_order_adjacent_vertex (const unsigned int n,
289  const unsigned int v) const
290 {
291  libmesh_assert_greater_equal (n, this->n_vertices());
292  libmesh_assert_less (n, this->n_nodes());
293 
294  switch (n)
295  {
296  case 5:
297  case 6:
298  case 7:
299  case 8:
300  case 9:
301  case 10:
302  case 11:
303  case 12:
304  {
305  libmesh_assert_less (v, 2);
306 
307  // This is the analog of the static, const arrays
308  // {Hex,Prism,Tet10}::_second_order_adjacent_vertices
309  // defined in the respective source files... possibly treat
310  // this similarly once the Pyramid13 has been added?
311  unsigned short node_list[8][2] =
312  {
313  {0,1},
314  {1,2},
315  {2,3},
316  {0,3},
317  {0,4},
318  {1,4},
319  {2,4},
320  {3,4}
321  };
322 
323  return node_list[n-5][v];
324  }
325 
326  // mid-face node on bottom
327  case 13:
328  {
329  libmesh_assert_less (v, 4);
330 
331  // The vertex nodes surrounding node 13 are 0, 1, 2, and 3.
332  // Thus, the v'th node is simply = v.
333  return cast_int<unsigned short>(v);
334  }
335 
336  default:
337  libmesh_error_msg("Invalid n = " << n);
338 
339  }
340 
341  libmesh_error_msg("We'll never get here!");
342  return static_cast<unsigned short int>(-1);
343 }
344 
345 
346 
347 Real Pyramid14::volume () const
348 {
349  // Make copies of our points. It makes the subsequent calculations a bit
350  // shorter and avoids dereferencing the same pointer multiple times.
351  Point
352  x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4), x5 = point(5), x6 = point(6),
353  x7 = point(7), x8 = point(8), x9 = point(9), x10 = point(10), x11 = point(11), x12 = point(12), x13 = point(13);
354 
355  // dx/dxi and dx/deta have 15 components while dx/dzeta has 20.
356  // These are copied directly from the output of a Python script.
357  Point dx_dxi[15] =
358  {
359  x6/2 - x8/2,
360  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - 3*x6/2 + 3*x8/2 - x9,
361  -x0/2 + x1/2 - 2*x10 - 2*x11 + 2*x12 + x2/2 - x3/2 + 3*x6/2 - 3*x8/2 + 2*x9,
362  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9,
363  x0/4 - x1/4 + x2/4 - x3/4,
364  -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9,
365  x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9,
366  -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2,
367  x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2,
368  -2*x13 + x6 + x8,
369  4*x13 - 2*x6 - 2*x8,
370  -2*x13 + x6 + x8,
371  -x0/2 - x1/2 + x2/2 + x3/2 + x5 - x7,
372  x0/2 + x1/2 - x2/2 - x3/2 - x5 + x7,
373  x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
374  };
375 
376  // dx/dxi and dx/deta have 15 components while dx/dzeta has 20.
377  // These are copied directly from the output of a Python script.
378  Point dx_deta[15] =
379  {
380  -x5/2 + x7/2,
381  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + 3*x5/2 - 3*x7/2 - x9,
382  -x0/2 - x1/2 + 2*x10 - 2*x11 - 2*x12 + x2/2 + x3/2 - 3*x5/2 + 3*x7/2 + 2*x9,
383  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
384  -2*x13 + x5 + x7,
385  4*x13 - 2*x5 - 2*x7,
386  -2*x13 + x5 + x7,
387  x0/4 - x1/4 + x2/4 - x3/4,
388  -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9,
389  x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9,
390  -x0/2 + x1/2 + x2/2 - x3/2 - x6 + x8,
391  x0/2 - x1/2 - x2/2 + x3/2 + x6 - x8,
392  -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
393  x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2,
394  x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
395  };
396 
397  // dx/dxi and dx/deta have 15 components while dx/dzeta has 20.
398  // These are copied directly from the output of a Python script.
399  Point dx_dzeta[20] =
400  {
401  -x0/4 - x1/4 + x10 + x11 + x12 - 2*x13 - x2/4 - x3/4 - x4 + x9,
402  5*x0/4 + 5*x1/4 - 5*x10 - 5*x11 - 5*x12 + 8*x13 + 5*x2/4 + 5*x3/4 + 7*x4 - 5*x9,
403  -9*x0/4 - 9*x1/4 + 9*x10 + 9*x11 + 9*x12 - 12*x13 - 9*x2/4 - 9*x3/4 - 15*x4 + 9*x9,
404  7*x0/4 + 7*x1/4 - 7*x10 - 7*x11 - 7*x12 + 8*x13 + 7*x2/4 + 7*x3/4 + 13*x4 - 7*x9,
405  -x0/2 - x1/2 + 2*x10 + 2*x11 + 2*x12 - 2*x13 - x2/2 - x3/2 - 4*x4 + 2*x9,
406  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
407  -3*x0/4 - 3*x1/4 + 3*x10 - 3*x11 - 3*x12 + 3*x2/4 + 3*x3/4 - 3*x5/2 + 3*x7/2 + 3*x9,
408  3*x0/4 + 3*x1/4 - 3*x10 + 3*x11 + 3*x12 - 3*x2/4 - 3*x3/4 + 3*x5/2 - 3*x7/2 - 3*x9,
409  -x0/4 - x1/4 + x10 - x11 - x12 + x2/4 + x3/4 - x5/2 + x7/2 + x9,
410  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9,
411  -3*x0/4 + 3*x1/4 - 3*x10 - 3*x11 + 3*x12 + 3*x2/4 - 3*x3/4 + 3*x6/2 - 3*x8/2 + 3*x9,
412  3*x0/4 - 3*x1/4 + 3*x10 + 3*x11 - 3*x12 - 3*x2/4 + 3*x3/4 - 3*x6/2 + 3*x8/2 - 3*x9,
413  -x0/4 + x1/4 - x10 - x11 + x12 + x2/4 - x3/4 + x6/2 - x8/2 + x9,
414  -x0/4 + x1/4 - x10 + x11 - x12 - x2/4 + x3/4 + x9,
415  x0/4 - x1/4 + x10 - x11 + x12 + x2/4 - x3/4 - x9,
416  -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2,
417  x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2,
418  -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
419  x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2,
420  x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
421  };
422 
423  // The (xi, eta, zeta) exponents for each of the dx_dxi terms
424  static const int dx_dxi_exponents[15][3] =
425  {
426  {0, 0, 0},
427  {0, 0, 1},
428  {0, 0, 2},
429  {0, 0, 3},
430  {0, 1, 0},
431  {0, 1, 1},
432  {0, 1, 2},
433  {0, 2, 0},
434  {0, 2, 1},
435  {1, 0, 0},
436  {1, 0, 1},
437  {1, 0, 2},
438  {1, 1, 0},
439  {1, 1, 1},
440  {1, 2, 0}
441  };
442 
443  // The (xi, eta, zeta) exponents for each of the dx_deta terms
444  static const int dx_deta_exponents[15][3] =
445  {
446  {0, 0, 0},
447  {0, 0, 1},
448  {0, 0, 2},
449  {0, 0, 3},
450  {0, 1, 0},
451  {0, 1, 1},
452  {0, 1, 2},
453  {1, 0, 0},
454  {1, 0, 1},
455  {1, 0, 2},
456  {1, 1, 0},
457  {1, 1, 1},
458  {2, 0, 0},
459  {2, 0, 1},
460  {2, 1, 0}
461  };
462 
463  // The (xi, eta, zeta) exponents for each of the dx_dzeta terms
464  static const int dx_dzeta_exponents[20][3] =
465  {
466  {0, 0, 0},
467  {0, 0, 1},
468  {0, 0, 2},
469  {0, 0, 3},
470  {0, 0, 4},
471  {0, 1, 0},
472  {0, 1, 1},
473  {0, 1, 2},
474  {0, 1, 3},
475  {1, 0, 0},
476  {1, 0, 1},
477  {1, 0, 2},
478  {1, 0, 3},
479  {1, 1, 0},
480  {1, 1, 1},
481  {1, 2, 0},
482  {1, 2, 1},
483  {2, 1, 0},
484  {2, 1, 1},
485  {2, 2, 0},
486  };
487 
488  // Number of points in the quadrature rule
489  const int N = 27;
490 
491  // Parameters of the quadrature rule
492  static const Real
493  // Parameters used for (xi, eta) quadrature points.
494  a1 = -7.1805574131988893873307823958101e-01L,
495  a2 = -5.0580870785392503961340276902425e-01L,
496  a3 = -2.2850430565396735359574351631314e-01L,
497  // Parameters used for zeta quadrature points.
498  b1 = 7.2994024073149732155837979012003e-02L,
499  b2 = 3.4700376603835188472176354340395e-01L,
500  b3 = 7.0500220988849838312239847758405e-01L,
501  // There are 9 unique weight values since there are three
502  // for each of the three unique zeta values.
503  w1 = 4.8498876871878584357513834016440e-02L,
504  w2 = 4.5137737425884574692441981593901e-02L,
505  w3 = 9.2440441384508327195915094925393e-03L,
506  w4 = 7.7598202995005734972022134426305e-02L,
507  w5 = 7.2220379881415319507907170550242e-02L,
508  w6 = 1.4790470621521332351346415188063e-02L,
509  w7 = 1.2415712479200917595523541508209e-01L,
510  w8 = 1.1555260781026451121265147288039e-01L,
511  w9 = 2.3664752994434131762154264300901e-02L;
512 
513  // The points and weights of the 3x3x3 quadrature rule
514  static const Real xi[N][3] =
515  {// ^0 ^1 ^2
516  { 1., a1, a1*a1},
517  { 1., a2, a2*a2},
518  { 1., a3, a3*a3},
519  { 1., a1, a1*a1},
520  { 1., a2, a2*a2},
521  { 1., a3, a3*a3},
522  { 1., a1, a1*a1},
523  { 1., a2, a2*a2},
524  { 1., a3, a3*a3},
525  { 1., 0., 0. },
526  { 1., 0., 0. },
527  { 1., 0., 0. },
528  { 1., 0., 0. },
529  { 1., 0., 0. },
530  { 1., 0., 0. },
531  { 1., 0., 0. },
532  { 1., 0., 0. },
533  { 1., 0., 0. },
534  { 1., -a1, a1*a1},
535  { 1., -a2, a2*a2},
536  { 1., -a3, a3*a3},
537  { 1., -a1, a1*a1},
538  { 1., -a2, a2*a2},
539  { 1., -a3, a3*a3},
540  { 1., -a1, a1*a1},
541  { 1., -a2, a2*a2},
542  { 1., -a3, a3*a3}
543  };
544 
545  static const Real eta[N][3] =
546  {// ^0 ^1 ^2
547  { 1., a1, a1*a1},
548  { 1., a2, a2*a2},
549  { 1., a3, a3*a3},
550  { 1., 0., 0. },
551  { 1., 0., 0. },
552  { 1., 0., 0. },
553  { 1., -a1, a1*a1},
554  { 1., -a2, a2*a2},
555  { 1., -a3, a3*a3},
556  { 1., a1, a1*a1},
557  { 1., a2, a2*a2},
558  { 1., a3, a3*a3},
559  { 1., 0., 0. },
560  { 1., 0., 0. },
561  { 1., 0., 0. },
562  { 1., -a1, a1*a1},
563  { 1., -a2, a2*a2},
564  { 1., -a3, a3*a3},
565  { 1., a1, a1*a1},
566  { 1., a2, a2*a2},
567  { 1., a3, a3*a3},
568  { 1., 0., 0. },
569  { 1., 0., 0. },
570  { 1., 0., 0. },
571  { 1., -a1, a1*a1},
572  { 1., -a2, a2*a2},
573  { 1., -a3, a3*a3}
574  };
575 
576  static const Real zeta[N][5] =
577  {// ^0 ^1 ^2 ^3 ^4
578  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
579  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
580  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
581  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
582  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
583  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
584  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
585  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
586  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
587  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
588  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
589  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
590  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
591  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
592  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
593  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
594  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
595  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
596  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
597  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
598  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
599  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
600  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
601  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
602  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
603  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
604  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}
605  };
606 
607  static const Real w[N] = {w1, w2, w3, w4, w5, w6, // 0-5
608  w1, w2, w3, w4, w5, w6, // 6-11
609  w7, w8, w9, w4, w5, w6, // 12-17
610  w1, w2, w3, w4, w5, w6, // 18-23
611  w1, w2, w3}; // 24-26
612 
613  Real vol = 0.;
614  for (int q=0; q<N; ++q)
615  {
616  // Compute denominators for the current q.
617  Real
618  den2 = (1. - zeta[q][1])*(1. - zeta[q][1]),
619  den3 = den2*(1. - zeta[q][1]);
620 
621  // Compute dx/dxi and dx/deta at the current q.
622  Point dx_dxi_q, dx_deta_q;
623  for (int c=0; c<15; ++c)
624  {
625  dx_dxi_q +=
626  xi[q][dx_dxi_exponents[c][0]]*
627  eta[q][dx_dxi_exponents[c][1]]*
628  zeta[q][dx_dxi_exponents[c][2]]*dx_dxi[c];
629 
630  dx_deta_q +=
631  xi[q][dx_deta_exponents[c][0]]*
632  eta[q][dx_deta_exponents[c][1]]*
633  zeta[q][dx_deta_exponents[c][2]]*dx_deta[c];
634  }
635 
636  // Compute dx/dzeta at the current q.
637  Point dx_dzeta_q;
638  for (int c=0; c<20; ++c)
639  {
640  dx_dzeta_q +=
641  xi[q][dx_dzeta_exponents[c][0]]*
642  eta[q][dx_dzeta_exponents[c][1]]*
643  zeta[q][dx_dzeta_exponents[c][2]]*dx_dzeta[c];
644  }
645 
646  // Scale everything appropriately
647  dx_dxi_q /= den2;
648  dx_deta_q /= den2;
649  dx_dzeta_q /= den3;
650 
651  // Compute scalar triple product, multiply by weight, and accumulate volume.
652  vol += w[q] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q);
653  }
654 
655  return vol;
656 }
657 
658 } // namespace libMesh
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:184
unsigned short int side
Definition: xdr_io.C:49
const class libmesh_nullptr_t libmesh_nullptr
The libMesh namespace provides an interface to certain functionality in the library.
libmesh_assert(j)
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1051
const dof_id_type n_nodes
Definition: tecplot_io.C:67
virtual dof_id_type key() const
Definition: elem.C:503
PetscErrorCode Vec Mat libmesh_dbg_var(j)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
uint8_t dof_id_type
Definition: id_types.h:64