libMesh
cell_pyramid14.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 // Local includes
20 #include "libmesh/cell_pyramid14.h"
21 #include "libmesh/edge_edge3.h"
22 #include "libmesh/face_tri6.h"
23 #include "libmesh/face_quad9.h"
24 #include "libmesh/enum_io_package.h"
25 #include "libmesh/enum_order.h"
26 
27 namespace libMesh
28 {
29 
30 
31 
32 
33 // ------------------------------------------------------------
34 // Pyramid14 class static member initializations
35 const int Pyramid14::num_nodes;
38 
40  {
41  {0, 1, 4, 5, 10, 9, 99, 99, 99}, // Side 0 (front)
42  {1, 2, 4, 6, 11, 10, 99, 99, 99}, // Side 1 (right)
43  {2, 3, 4, 7, 12, 11, 99, 99, 99}, // Side 2 (back)
44  {3, 0, 4, 8, 9, 12, 99, 99, 99}, // Side 3 (left)
45  {0, 3, 2, 1, 8, 7, 6, 5, 13} // Side 4 (base)
46  };
47 
49  {
50  {0, 1, 5}, // Edge 0
51  {1, 2, 6}, // Edge 1
52  {2, 3, 7}, // Edge 2
53  {0, 3, 8}, // Edge 3
54  {0, 4, 9}, // Edge 4
55  {1, 4, 10}, // Edge 5
56  {2, 4, 11}, // Edge 6
57  {3, 4, 12} // Edge 7
58  };
59 
60 // ------------------------------------------------------------
61 // Pyramid14 class member functions
62 
63 bool Pyramid14::is_vertex(const unsigned int i) const
64 {
65  if (i < 5)
66  return true;
67  return false;
68 }
69 
70 
71 
72 bool Pyramid14::is_edge(const unsigned int i) const
73 {
74  if (i < 5)
75  return false;
76  if (i == 13)
77  return false;
78  return true;
79 }
80 
81 
82 
83 bool Pyramid14::is_face(const unsigned int i) const
84 {
85  if (i == 13)
86  return true;
87  return false;
88 }
89 
90 
91 
92 bool Pyramid14::is_node_on_side(const unsigned int n,
93  const unsigned int s) const
94 {
95  libmesh_assert_less (s, n_sides());
96  return std::find(std::begin(side_nodes_map[s]),
97  std::end(side_nodes_map[s]),
98  n) != std::end(side_nodes_map[s]);
99 }
100 
101 std::vector<unsigned>
102 Pyramid14::nodes_on_side(const unsigned int s) const
103 {
104  libmesh_assert_less(s, n_sides());
105  auto trim = (s == 4) ? 0 : 3;
106  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s]) - trim};
107 }
108 
109 std::vector<unsigned>
110 Pyramid14::nodes_on_edge(const unsigned int e) const
111 {
112  libmesh_assert_less(e, n_edges());
113  return {std::begin(edge_nodes_map[e]), std::end(edge_nodes_map[e])};
114 }
115 
116 bool Pyramid14::is_node_on_edge(const unsigned int n,
117  const unsigned int e) const
118 {
119  libmesh_assert_less (e, n_edges());
120  return std::find(std::begin(edge_nodes_map[e]),
121  std::end(edge_nodes_map[e]),
122  n) != std::end(edge_nodes_map[e]);
123 }
124 
126 {
127  // TODO: If the base is a parallelogram and all the triangular faces are planar,
128  // the map should be linear, but I need to test this theory...
129  return false;
130 }
131 
132 
133 
135 {
136  return SECOND;
137 }
138 
139 
140 
141 dof_id_type Pyramid14::key (const unsigned int s) const
142 {
143  libmesh_assert_less (s, this->n_sides());
144 
145  switch (s)
146  {
147  case 0: // triangular face 1
148  case 1: // triangular face 2
149  case 2: // triangular face 3
150  case 3: // triangular face 4
151  return Pyramid::key(s);
152 
153  case 4: // the quad face at z=0
154  return this->compute_key (this->node_id(13));
155 
156  default:
157  libmesh_error_msg("Invalid side s = " << s);
158  }
159 }
160 
161 
162 
163 unsigned int Pyramid14::local_side_node(unsigned int side,
164  unsigned int side_node) const
165 {
166  libmesh_assert_less (side, this->n_sides());
167 
168  // Never more than 9 nodes per side.
169  libmesh_assert_less(side_node, Pyramid14::nodes_per_side);
170 
171  // Some sides have 6 nodes.
172  libmesh_assert(side == 4 || side_node < 6);
173 
174  return Pyramid14::side_nodes_map[side][side_node];
175 }
176 
177 
178 
179 unsigned int Pyramid14::local_edge_node(unsigned int edge,
180  unsigned int edge_node) const
181 {
182  libmesh_assert_less(edge, this->n_edges());
183  libmesh_assert_less(edge_node, Pyramid14::nodes_per_edge);
184 
185  return Pyramid14::edge_nodes_map[edge][edge_node];
186 }
187 
188 
189 
190 std::unique_ptr<Elem> Pyramid14::build_side_ptr (const unsigned int i)
191 {
192  libmesh_assert_less (i, this->n_sides());
193 
194  std::unique_ptr<Elem> face;
195 
196  switch (i)
197  {
198  case 0: // triangular face 1
199  case 1: // triangular face 2
200  case 2: // triangular face 3
201  case 3: // triangular face 4
202  {
203  face = std::make_unique<Tri6>();
204  break;
205  }
206  case 4: // the quad face at z=0
207  {
208  face = std::make_unique<Quad9>();
209  break;
210  }
211  default:
212  libmesh_error_msg("Invalid side i = " << i);
213  }
214 
215  // Set the nodes
216  for (auto n : face->node_index_range())
217  face->set_node(n, this->node_ptr(Pyramid14::side_nodes_map[i][n]));
218 
219  face->set_interior_parent(this);
220  face->inherit_data_from(*this);
221 
222  return face;
223 }
224 
225 
226 
227 void Pyramid14::build_side_ptr (std::unique_ptr<Elem> & side,
228  const unsigned int i)
229 {
230  libmesh_assert_less (i, this->n_sides());
231 
232  switch (i)
233  {
234  case 0: // triangular face 1
235  case 1: // triangular face 2
236  case 2: // triangular face 3
237  case 3: // triangular face 4
238  {
239  if (!side.get() || side->type() != TRI6)
240  {
241  side = this->build_side_ptr(i);
242  return;
243  }
244  break;
245  }
246  case 4: // the quad face at z=0
247  {
248  if (!side.get() || side->type() != QUAD9)
249  {
250  side = this->build_side_ptr(i);
251  return;
252  }
253  break;
254  }
255  default:
256  libmesh_error_msg("Invalid side i = " << i);
257  }
258 
259  side->inherit_data_from(*this);
260 
261  // Set the nodes
262  for (auto n : side->node_index_range())
263  side->set_node(n, this->node_ptr(Pyramid14::side_nodes_map[i][n]));
264 }
265 
266 
267 
268 std::unique_ptr<Elem> Pyramid14::build_edge_ptr (const unsigned int i)
269 {
270  return this->simple_build_edge_ptr<Edge3,Pyramid14>(i);
271 }
272 
273 
274 
275 void Pyramid14::build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i)
276 {
277  this->simple_build_edge_ptr<Pyramid14>(edge, i, EDGE3);
278 }
279 
280 
281 
282 void Pyramid14::connectivity(const unsigned int libmesh_dbg_var(sc),
283  const IOPackage iop,
284  std::vector<dof_id_type> & /*conn*/) const
285 {
287  libmesh_assert_less (sc, this->n_sub_elem());
288  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
289 
290  switch (iop)
291  {
292  case TECPLOT:
293  {
294  // TODO
295  libmesh_not_implemented();
296  }
297 
298  case VTK:
299  {
300  // TODO
301  libmesh_not_implemented();
302  }
303 
304  default:
305  libmesh_error_msg("Unsupported IO package " << iop);
306  }
307 }
308 
309 
310 
311 unsigned int Pyramid14::n_second_order_adjacent_vertices (const unsigned int n) const
312 {
313  switch (n)
314  {
315  case 5:
316  case 6:
317  case 7:
318  case 8:
319  case 9:
320  case 10:
321  case 11:
322  case 12:
323  return 2;
324 
325  case 13:
326  return 4;
327 
328  default:
329  libmesh_error_msg("Invalid node n = " << n);
330  }
331 }
332 
333 
334 unsigned short int Pyramid14::second_order_adjacent_vertex (const unsigned int n,
335  const unsigned int v) const
336 {
337  libmesh_assert_greater_equal (n, this->n_vertices());
338  libmesh_assert_less (n, this->n_nodes());
339 
340  switch (n)
341  {
342  case 5:
343  case 6:
344  case 7:
345  case 8:
346  case 9:
347  case 10:
348  case 11:
349  case 12:
350  {
351  libmesh_assert_less (v, 2);
352 
353  // This is the analog of the static, const arrays
354  // {Hex,Prism,Tet10}::_second_order_adjacent_vertices
355  // defined in the respective source files... possibly treat
356  // this similarly once the Pyramid13 has been added?
357  unsigned short node_list[8][2] =
358  {
359  {0,1},
360  {1,2},
361  {2,3},
362  {0,3},
363  {0,4},
364  {1,4},
365  {2,4},
366  {3,4}
367  };
368 
369  return node_list[n-5][v];
370  }
371 
372  // mid-face node on bottom
373  case 13:
374  {
375  libmesh_assert_less (v, 4);
376 
377  // The vertex nodes surrounding node 13 are 0, 1, 2, and 3.
378  // Thus, the v'th node is simply = v.
379  return cast_int<unsigned short>(v);
380  }
381 
382  default:
383  libmesh_error_msg("Invalid n = " << n);
384 
385  }
386 }
387 
388 
389 
391 {
392  // This specialization is good for Lagrange mappings only
393  if (this->mapping_type() != LAGRANGE_MAP)
394  return this->Elem::volume();
395 
396  // Make copies of our points. It makes the subsequent calculations a bit
397  // shorter and avoids dereferencing the same pointer multiple times.
398  Point
399  x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4), x5 = point(5), x6 = point(6),
400  x7 = point(7), x8 = point(8), x9 = point(9), x10 = point(10), x11 = point(11), x12 = point(12), x13 = point(13);
401 
402  // dx/dxi and dx/deta have 15 components while dx/dzeta has 20.
403  // These are copied directly from the output of a Python script.
404  Point dx_dxi[15] =
405  {
406  x6/2 - x8/2,
407  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - 3*x6/2 + 3*x8/2 - x9,
408  -x0/2 + x1/2 - 2*x10 - 2*x11 + 2*x12 + x2/2 - x3/2 + 3*x6/2 - 3*x8/2 + 2*x9,
409  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9,
410  x0/4 - x1/4 + x2/4 - x3/4,
411  -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9,
412  x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9,
413  -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2,
414  x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2,
415  -2*x13 + x6 + x8,
416  4*x13 - 2*x6 - 2*x8,
417  -2*x13 + x6 + x8,
418  -x0/2 - x1/2 + x2/2 + x3/2 + x5 - x7,
419  x0/2 + x1/2 - x2/2 - x3/2 - x5 + x7,
420  x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
421  };
422 
423  // dx/dxi and dx/deta have 15 components while dx/dzeta has 20.
424  // These are copied directly from the output of a Python script.
425  Point dx_deta[15] =
426  {
427  -x5/2 + x7/2,
428  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + 3*x5/2 - 3*x7/2 - x9,
429  -x0/2 - x1/2 + 2*x10 - 2*x11 - 2*x12 + x2/2 + x3/2 - 3*x5/2 + 3*x7/2 + 2*x9,
430  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
431  -2*x13 + x5 + x7,
432  4*x13 - 2*x5 - 2*x7,
433  -2*x13 + x5 + x7,
434  x0/4 - x1/4 + x2/4 - x3/4,
435  -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9,
436  x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9,
437  -x0/2 + x1/2 + x2/2 - x3/2 - x6 + x8,
438  x0/2 - x1/2 - x2/2 + x3/2 + x6 - x8,
439  -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
440  x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2,
441  x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
442  };
443 
444  // dx/dxi and dx/deta have 15 components while dx/dzeta has 20.
445  // These are copied directly from the output of a Python script.
446  Point dx_dzeta[20] =
447  {
448  -x0/4 - x1/4 + x10 + x11 + x12 - 2*x13 - x2/4 - x3/4 - x4 + x9,
449  5*x0/4 + 5*x1/4 - 5*x10 - 5*x11 - 5*x12 + 8*x13 + 5*x2/4 + 5*x3/4 + 7*x4 - 5*x9,
450  -9*x0/4 - 9*x1/4 + 9*x10 + 9*x11 + 9*x12 - 12*x13 - 9*x2/4 - 9*x3/4 - 15*x4 + 9*x9,
451  7*x0/4 + 7*x1/4 - 7*x10 - 7*x11 - 7*x12 + 8*x13 + 7*x2/4 + 7*x3/4 + 13*x4 - 7*x9,
452  -x0/2 - x1/2 + 2*x10 + 2*x11 + 2*x12 - 2*x13 - x2/2 - x3/2 - 4*x4 + 2*x9,
453  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
454  -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,
455  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,
456  -x0/4 - x1/4 + x10 - x11 - x12 + x2/4 + x3/4 - x5/2 + x7/2 + x9,
457  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9,
458  -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,
459  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,
460  -x0/4 + x1/4 - x10 - x11 + x12 + x2/4 - x3/4 + x6/2 - x8/2 + x9,
461  -x0/4 + x1/4 - x10 + x11 - x12 - x2/4 + x3/4 + x9,
462  x0/4 - x1/4 + x10 - x11 + x12 + x2/4 - x3/4 - x9,
463  -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2,
464  x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2,
465  -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
466  x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2,
467  x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
468  };
469 
470  // The (xi, eta, zeta) exponents for each of the dx_dxi terms
471  static const int dx_dxi_exponents[15][3] =
472  {
473  {0, 0, 0},
474  {0, 0, 1},
475  {0, 0, 2},
476  {0, 0, 3},
477  {0, 1, 0},
478  {0, 1, 1},
479  {0, 1, 2},
480  {0, 2, 0},
481  {0, 2, 1},
482  {1, 0, 0},
483  {1, 0, 1},
484  {1, 0, 2},
485  {1, 1, 0},
486  {1, 1, 1},
487  {1, 2, 0}
488  };
489 
490  // The (xi, eta, zeta) exponents for each of the dx_deta terms
491  static const int dx_deta_exponents[15][3] =
492  {
493  {0, 0, 0},
494  {0, 0, 1},
495  {0, 0, 2},
496  {0, 0, 3},
497  {0, 1, 0},
498  {0, 1, 1},
499  {0, 1, 2},
500  {1, 0, 0},
501  {1, 0, 1},
502  {1, 0, 2},
503  {1, 1, 0},
504  {1, 1, 1},
505  {2, 0, 0},
506  {2, 0, 1},
507  {2, 1, 0}
508  };
509 
510  // The (xi, eta, zeta) exponents for each of the dx_dzeta terms
511  static const int dx_dzeta_exponents[20][3] =
512  {
513  {0, 0, 0},
514  {0, 0, 1},
515  {0, 0, 2},
516  {0, 0, 3},
517  {0, 0, 4},
518  {0, 1, 0},
519  {0, 1, 1},
520  {0, 1, 2},
521  {0, 1, 3},
522  {1, 0, 0},
523  {1, 0, 1},
524  {1, 0, 2},
525  {1, 0, 3},
526  {1, 1, 0},
527  {1, 1, 1},
528  {1, 2, 0},
529  {1, 2, 1},
530  {2, 1, 0},
531  {2, 1, 1},
532  {2, 2, 0},
533  };
534 
535  // Number of points in the quadrature rule
536  const int N = 27;
537 
538  // Parameters of the quadrature rule
539  static const Real
540  // Parameters used for (xi, eta) quadrature points.
541  a1 = -7.1805574131988893873307823958101e-01_R,
542  a2 = -5.0580870785392503961340276902425e-01_R,
543  a3 = -2.2850430565396735359574351631314e-01_R,
544  // Parameters used for zeta quadrature points.
545  b1 = 7.2994024073149732155837979012003e-02_R,
546  b2 = 3.4700376603835188472176354340395e-01_R,
547  b3 = 7.0500220988849838312239847758405e-01_R,
548  // There are 9 unique weight values since there are three
549  // for each of the three unique zeta values.
550  w1 = 4.8498876871878584357513834016440e-02_R,
551  w2 = 4.5137737425884574692441981593901e-02_R,
552  w3 = 9.2440441384508327195915094925393e-03_R,
553  w4 = 7.7598202995005734972022134426305e-02_R,
554  w5 = 7.2220379881415319507907170550242e-02_R,
555  w6 = 1.4790470621521332351346415188063e-02_R,
556  w7 = 1.2415712479200917595523541508209e-01_R,
557  w8 = 1.1555260781026451121265147288039e-01_R,
558  w9 = 2.3664752994434131762154264300901e-02_R;
559 
560  // The points and weights of the 3x3x3 quadrature rule
561  static const Real xi[N][3] =
562  {// ^0 ^1 ^2
563  { 1., a1, a1*a1},
564  { 1., a2, a2*a2},
565  { 1., a3, a3*a3},
566  { 1., a1, a1*a1},
567  { 1., a2, a2*a2},
568  { 1., a3, a3*a3},
569  { 1., a1, a1*a1},
570  { 1., a2, a2*a2},
571  { 1., a3, a3*a3},
572  { 1., 0., 0. },
573  { 1., 0., 0. },
574  { 1., 0., 0. },
575  { 1., 0., 0. },
576  { 1., 0., 0. },
577  { 1., 0., 0. },
578  { 1., 0., 0. },
579  { 1., 0., 0. },
580  { 1., 0., 0. },
581  { 1., -a1, a1*a1},
582  { 1., -a2, a2*a2},
583  { 1., -a3, a3*a3},
584  { 1., -a1, a1*a1},
585  { 1., -a2, a2*a2},
586  { 1., -a3, a3*a3},
587  { 1., -a1, a1*a1},
588  { 1., -a2, a2*a2},
589  { 1., -a3, a3*a3}
590  };
591 
592  static const Real eta[N][3] =
593  {// ^0 ^1 ^2
594  { 1., a1, a1*a1},
595  { 1., a2, a2*a2},
596  { 1., a3, a3*a3},
597  { 1., 0., 0. },
598  { 1., 0., 0. },
599  { 1., 0., 0. },
600  { 1., -a1, a1*a1},
601  { 1., -a2, a2*a2},
602  { 1., -a3, a3*a3},
603  { 1., a1, a1*a1},
604  { 1., a2, a2*a2},
605  { 1., a3, a3*a3},
606  { 1., 0., 0. },
607  { 1., 0., 0. },
608  { 1., 0., 0. },
609  { 1., -a1, a1*a1},
610  { 1., -a2, a2*a2},
611  { 1., -a3, a3*a3},
612  { 1., a1, a1*a1},
613  { 1., a2, a2*a2},
614  { 1., a3, a3*a3},
615  { 1., 0., 0. },
616  { 1., 0., 0. },
617  { 1., 0., 0. },
618  { 1., -a1, a1*a1},
619  { 1., -a2, a2*a2},
620  { 1., -a3, a3*a3}
621  };
622 
623  static const Real zeta[N][5] =
624  {// ^0 ^1 ^2 ^3 ^4
625  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
626  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
627  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
628  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
629  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
630  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
631  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
632  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
633  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
634  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
635  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
636  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
637  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
638  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
639  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
640  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
641  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
642  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
643  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
644  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
645  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
646  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
647  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
648  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
649  { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
650  { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
651  { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}
652  };
653 
654  static const Real w[N] = {w1, w2, w3, w4, w5, w6, // 0-5
655  w1, w2, w3, w4, w5, w6, // 6-11
656  w7, w8, w9, w4, w5, w6, // 12-17
657  w1, w2, w3, w4, w5, w6, // 18-23
658  w1, w2, w3}; // 24-26
659 
660  Real vol = 0.;
661  for (int q=0; q<N; ++q)
662  {
663  // Compute denominators for the current q.
664  Real
665  den2 = (1. - zeta[q][1])*(1. - zeta[q][1]),
666  den3 = den2*(1. - zeta[q][1]);
667 
668  // Compute dx/dxi and dx/deta at the current q.
669  Point dx_dxi_q, dx_deta_q;
670  for (int c=0; c<15; ++c)
671  {
672  dx_dxi_q +=
673  xi[q][dx_dxi_exponents[c][0]]*
674  eta[q][dx_dxi_exponents[c][1]]*
675  zeta[q][dx_dxi_exponents[c][2]]*dx_dxi[c];
676 
677  dx_deta_q +=
678  xi[q][dx_deta_exponents[c][0]]*
679  eta[q][dx_deta_exponents[c][1]]*
680  zeta[q][dx_deta_exponents[c][2]]*dx_deta[c];
681  }
682 
683  // Compute dx/dzeta at the current q.
684  Point dx_dzeta_q;
685  for (int c=0; c<20; ++c)
686  {
687  dx_dzeta_q +=
688  xi[q][dx_dzeta_exponents[c][0]]*
689  eta[q][dx_dzeta_exponents[c][1]]*
690  zeta[q][dx_dzeta_exponents[c][2]]*dx_dzeta[c];
691  }
692 
693  // Scale everything appropriately
694  dx_dxi_q /= den2;
695  dx_deta_q /= den2;
696  dx_dzeta_q /= den3;
697 
698  // Compute scalar triple product, multiply by weight, and accumulate volume.
699  vol += w[q] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q);
700  }
701 
702  return vol;
703 }
704 
705 
706 void Pyramid14::permute(unsigned int perm_num)
707 {
708  libmesh_assert_less (perm_num, 4);
709 
710  for (unsigned int i = 0; i != perm_num; ++i)
711  {
712  swap4nodes(0,1,2,3);
713  swap4nodes(5,6,7,8);
714  swap4nodes(9,10,11,12);
715  swap4neighbors(0,1,2,3);
716  }
717 }
718 
719 
720 void Pyramid14::flip(BoundaryInfo * boundary_info)
721 {
722  libmesh_assert(boundary_info);
723 
724  swap2nodes(0,1);
725  swap2nodes(2,3);
726  swap2nodes(6,8);
727  swap2nodes(9,10);
728  swap2nodes(11,12);
729  swap2neighbors(1,3);
730  swap2boundarysides(1,3,boundary_info);
731  swap2boundaryedges(1,3,boundary_info);
732  swap2boundaryedges(4,5,boundary_info);
733  swap2boundaryedges(6,7,boundary_info);
734 }
735 
736 
737 unsigned int Pyramid14::center_node_on_side(const unsigned short side) const
738 {
739  libmesh_assert_less (side, Pyramid14::num_sides);
740  return side == 4 ? 13 : invalid_uint;
741 }
742 
743 
744 ElemType Pyramid14::side_type (const unsigned int s) const
745 {
746  libmesh_assert_less (s, 5);
747  if (s < 4)
748  return TRI6;
749  return QUAD9;
750 }
751 
752 
753 } // namespace libMesh
ElemType
Defines an enum for geometric element types.
void swap2boundaryedges(unsigned short e1, unsigned short e2, BoundaryInfo *boundary_info) const
Swaps two edges in boundary_info, if it is non-null.
Definition: elem.C:3550
virtual unsigned short int second_order_adjacent_vertex(const unsigned int n, const unsigned int v) const override
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
Node ** _nodes
Pointers to the nodes we are connected to.
Definition: elem.h:2245
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
Definition: libmesh.h:310
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
static const int num_edges
Definition: cell_pyramid.h:78
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override
virtual bool has_affine_map() const override
virtual dof_id_type key() const
Definition: elem.C:753
virtual Order default_order() const override
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
void swap2boundarysides(unsigned short s1, unsigned short s2, BoundaryInfo *boundary_info) const
Swaps two sides in boundary_info, if it is non-null.
Definition: elem.C:3534
virtual unsigned int n_sides() const override
Definition: cell_pyramid.h:90
The libMesh namespace provides an interface to certain functionality in the library.
virtual unsigned int n_vertices() const override
Definition: cell_pyramid.h:95
virtual bool is_face(const unsigned int i) const override
virtual bool is_edge(const unsigned int i) const override
virtual unsigned int n_sub_elem() const override
FIXME: we don&#39;t yet have a refinement pattern for pyramids...
virtual void permute(unsigned int perm_num) override final
Permutes the element (by swapping node and neighbor pointers) according to the specified index...
static const int nodes_per_side
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1029
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
void swap4nodes(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four node_ptrs, "rotating" them.
Definition: elem.h:2143
ElemMappingType mapping_type() const
Definition: elem.h:3120
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:2092
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
static const int num_sides
Geometric constants for all Pyramids.
Definition: cell_pyramid.h:77
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE3 coincident with edge i.
static const int num_nodes
Geometric constants for Pyramid14.
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:2102
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned int n_nodes() const override
virtual bool is_vertex(const unsigned int i) const override
virtual Real volume() const override
Specialization for computing the volume of a Pyramid14.
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i) override
Builds a QUAD9 or TRI6 coincident with face i.
virtual unsigned int n_edges() const override
Definition: cell_pyramid.h:100
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
This maps the node of the edge to element node numbers.
virtual Real volume() const
Definition: elem.C:3429
void swap4neighbors(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four neighbor_ptrs, "rotating" them.
Definition: elem.h:2153
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
static dof_id_type compute_key(dof_id_type n0)
Definition: elem.h:3294
virtual unsigned int n_second_order_adjacent_vertices(const unsigned int n) const override
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
unsigned int center_node_on_side(const unsigned short side) const override final
ElemType side_type(const unsigned int s) const override final
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
dof_id_type node_id(const unsigned int i) const
Definition: elem.h:2475
const Point & point(const unsigned int i) const
Definition: elem.h:2453
static const int nodes_per_edge
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
uint8_t dof_id_type
Definition: id_types.h:67
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override