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