libMesh
cell_pyramid13.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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/side.h"
21 #include "libmesh/cell_pyramid13.h"
22 #include "libmesh/edge_edge3.h"
23 #include "libmesh/face_tri6.h"
24 #include "libmesh/face_quad8.h"
25 #include "libmesh/enum_io_package.h"
26 #include "libmesh/enum_order.h"
27 
28 namespace libMesh
29 {
30 
31 
32 
33 
34 // ------------------------------------------------------------
35 // Pyramid13 class static member initializations
36 const int Pyramid13::num_nodes;
37 const int Pyramid13::num_sides;
38 const int Pyramid13::num_edges;
39 const int Pyramid13::num_children;
42 
44  {
45  {0, 1, 4, 5, 10, 9, 99, 99}, // Side 0 (front)
46  {1, 2, 4, 6, 11, 10, 99, 99}, // Side 1 (right)
47  {2, 3, 4, 7, 12, 11, 99, 99}, // Side 2 (back)
48  {3, 0, 4, 8, 9, 12, 99, 99}, // Side 3 (left)
49  {0, 3, 2, 1, 8, 7, 6, 5} // Side 4 (base)
50  };
51 
53  {
54  {0, 1, 5}, // Edge 0
55  {1, 2, 6}, // Edge 1
56  {2, 3, 7}, // Edge 2
57  {0, 3, 8}, // Edge 3
58  {0, 4, 9}, // Edge 4
59  {1, 4, 10}, // Edge 5
60  {2, 4, 11}, // Edge 6
61  {3, 4, 12} // Edge 7
62  };
63 
64 // ------------------------------------------------------------
65 // Pyramid13 class member functions
66 
67 bool Pyramid13::is_vertex(const unsigned int i) const
68 {
69  if (i < 5)
70  return true;
71  return false;
72 }
73 
74 
75 
76 bool Pyramid13::is_edge(const unsigned int i) const
77 {
78  if (i < 5)
79  return false;
80  return true;
81 }
82 
83 
84 
85 bool Pyramid13::is_face(const unsigned int) const
86 {
87  return false;
88 }
89 
90 
91 
92 bool Pyramid13::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 Pyramid13::nodes_on_side(const unsigned int s) const
103 {
104  libmesh_assert_less(s, n_sides());
105  auto trim = (s == 4) ? 0 : 2;
106  return {std::begin(side_nodes_map[s]), std::end(side_nodes_map[s]) - trim};
107 }
108 
109 std::vector<unsigned>
110 Pyramid13::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 Pyramid13::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 
125 
126 
128 {
129  // TODO: If the base is a parallelogram and all the triangular faces are planar,
130  // the map should be linear, but I need to test this theory...
131  return false;
132 }
133 
134 
135 
137 {
138  return SECOND;
139 }
140 
141 
142 
143 unsigned int Pyramid13::local_side_node(unsigned int side,
144  unsigned int side_node) const
145 {
146  libmesh_assert_less (side, this->n_sides());
147 
148  // Never more than 8 nodes per side.
149  libmesh_assert_less(side_node, Pyramid13::nodes_per_side);
150 
151  // Some sides have 6 nodes.
152  libmesh_assert(side == 4 || side_node < 6);
153 
154  return Pyramid13::side_nodes_map[side][side_node];
155 }
156 
157 
158 
159 unsigned int Pyramid13::local_edge_node(unsigned int edge,
160  unsigned int edge_node) const
161 {
162  libmesh_assert_less(edge, this->n_edges());
163  libmesh_assert_less(edge_node, Pyramid13::nodes_per_edge);
164 
165  return Pyramid13::edge_nodes_map[edge][edge_node];
166 }
167 
168 
169 
170 std::unique_ptr<Elem> Pyramid13::build_side_ptr (const unsigned int i, bool proxy)
171 {
172  libmesh_assert_less (i, this->n_sides());
173 
174  std::unique_ptr<Elem> face;
175  if (proxy)
176  {
177 #ifdef LIBMESH_ENABLE_DEPRECATED
178  libmesh_deprecated();
179  switch (i)
180  {
181  case 0:
182  case 1:
183  case 2:
184  case 3:
185  {
186  face = std::make_unique<Side<Tri6,Pyramid13>>(this,i);
187  break;
188  }
189 
190  case 4:
191  {
192  face = std::make_unique<Side<Quad8,Pyramid13>>(this,i);
193  break;
194  }
195 
196  default:
197  libmesh_error_msg("Invalid side i = " << i);
198  }
199 #else
200  libmesh_error();
201 #endif // LIBMESH_ENABLE_DEPRECATED
202  }
203  else
204  {
205  switch (i)
206  {
207  case 0: // triangular face 1
208  case 1: // triangular face 2
209  case 2: // triangular face 3
210  case 3: // triangular face 4
211  {
212  face = std::make_unique<Tri6>();
213  break;
214  }
215  case 4: // the quad face at z=0
216  {
217  face = std::make_unique<Quad8>();
218  break;
219  }
220  default:
221  libmesh_error_msg("Invalid side i = " << i);
222  }
223 
224  // Set the nodes
225  for (auto n : face->node_index_range())
226  face->set_node(n) = this->node_ptr(Pyramid13::side_nodes_map[i][n]);
227  }
228 
229 #ifdef LIBMESH_ENABLE_DEPRECATED
230  if (!proxy) // proxy sides used to leave parent() set
231 #endif
232  face->set_parent(nullptr);
233  face->set_interior_parent(this);
234 
235  face->subdomain_id() = this->subdomain_id();
236  face->set_mapping_type(this->mapping_type());
237 #ifdef LIBMESH_ENABLE_AMR
238  face->set_p_level(this->p_level());
239 #endif
240 
241  return face;
242 }
243 
244 
245 
246 void Pyramid13::build_side_ptr (std::unique_ptr<Elem> & side,
247  const unsigned int i)
248 {
249  libmesh_assert_less (i, this->n_sides());
250 
251  switch (i)
252  {
253  case 0: // triangular face 1
254  case 1: // triangular face 2
255  case 2: // triangular face 3
256  case 3: // triangular face 4
257  {
258  if (!side.get() || side->type() != TRI6)
259  {
260  side = this->build_side_ptr(i, false);
261  return;
262  }
263  break;
264  }
265  case 4: // the quad face at z=0
266  {
267  if (!side.get() || side->type() != QUAD8)
268  {
269  side = this->build_side_ptr(i, false);
270  return;
271  }
272  break;
273  }
274  default:
275  libmesh_error_msg("Invalid side i = " << i);
276  }
277 
278  side->subdomain_id() = this->subdomain_id();
279  side->set_mapping_type(this->mapping_type());
280 
281  // Set the nodes
282  for (auto n : side->node_index_range())
283  side->set_node(n) = this->node_ptr(Pyramid13::side_nodes_map[i][n]);
284 }
285 
286 
287 
288 std::unique_ptr<Elem> Pyramid13::build_edge_ptr (const unsigned int i)
289 {
290  return this->simple_build_edge_ptr<Edge3,Pyramid13>(i);
291 }
292 
293 
294 
295 void Pyramid13::build_edge_ptr (std::unique_ptr<Elem> & edge, const unsigned int i)
296 {
297  this->simple_build_edge_ptr<Pyramid13>(edge, i, EDGE3);
298 }
299 
300 
301 
302 void Pyramid13::connectivity(const unsigned int libmesh_dbg_var(sc),
303  const IOPackage iop,
304  std::vector<dof_id_type> & /*conn*/) const
305 {
307  libmesh_assert_less (sc, this->n_sub_elem());
308  libmesh_assert_not_equal_to (iop, INVALID_IO_PACKAGE);
309 
310  switch (iop)
311  {
312  case TECPLOT:
313  {
314  // TODO
315  libmesh_not_implemented();
316  }
317 
318  case VTK:
319  {
320  // TODO
321  libmesh_not_implemented();
322  }
323 
324  default:
325  libmesh_error_msg("Unsupported IO package " << iop);
326  }
327 }
328 
329 
330 
331 unsigned int Pyramid13::n_second_order_adjacent_vertices (const unsigned int n) const
332 {
333  switch (n)
334  {
335  case 5:
336  case 6:
337  case 7:
338  case 8:
339  case 9:
340  case 10:
341  case 11:
342  case 12:
343  return 2;
344 
345  default:
346  libmesh_error_msg("Invalid node n = " << n);
347  }
348 }
349 
350 
351 unsigned short int Pyramid13::second_order_adjacent_vertex (const unsigned int n,
352  const unsigned int v) const
353 {
354  libmesh_assert_greater_equal (n, this->n_vertices());
355  libmesh_assert_less (n, this->n_nodes());
356 
357  switch (n)
358  {
359  case 5:
360  case 6:
361  case 7:
362  case 8:
363  case 9:
364  case 10:
365  case 11:
366  case 12:
367  {
368  libmesh_assert_less (v, 2);
369 
370  // This is the analog of the static, const arrays
371  // {Hex,Prism,Tet10}::_second_order_adjacent_vertices
372  // defined in the respective source files...
373  unsigned short node_list[8][2] =
374  {
375  {0,1},
376  {1,2},
377  {2,3},
378  {0,3},
379  {0,4},
380  {1,4},
381  {2,4},
382  {3,4}
383  };
384 
385  return node_list[n-5][v];
386  }
387 
388  default:
389  libmesh_error_msg("Invalid n = " << n);
390 
391  }
392 }
393 
394 
395 
397 {
398  // This specialization is good for Lagrange mappings only
399  if (this->mapping_type() != LAGRANGE_MAP)
400  return this->Elem::volume();
401 
402  // Make copies of our points. It makes the subsequent calculations a bit
403  // shorter and avoids dereferencing the same pointer multiple times.
404  Point
405  x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4), x5 = point(5), x6 = point(6),
406  x7 = point(7), x8 = point(8), x9 = point(9), x10 = point(10), x11 = point(11), x12 = point(12);
407 
408  // dx/dxi and dx/deta have 14 components while dx/dzeta has 19.
409  // These are copied directly from the output of a Python script.
410  Point dx_dxi[14] =
411  {
412  x6/2 - x8/2,
413  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - 3*x6/2 + 3*x8/2 - x9,
414  -x0/2 + x1/2 - 2*x10 - 2*x11 + 2*x12 + x2/2 - x3/2 + 3*x6/2 - 3*x8/2 + 2*x9,
415  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9,
416  x0/4 - x1/4 + x2/4 - x3/4,
417  -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9,
418  x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9,
419  -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2,
420  x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2,
421  x0/2 + x1/2 + x2/2 + x3/2 - x5 - x7,
422  -x0 - x1 - x2 - x3 + 2*x5 + 2*x7,
423  x0/2 + x1/2 + x2/2 + x3/2 - x5 - x7,
424  -x0/2 - x1/2 + x2/2 + x3/2 + x5 - x7,
425  x0/2 + x1/2 - x2/2 - x3/2 - x5 + x7
426  };
427 
428  // dx/dxi and dx/deta have 14 components while dx/dzeta has 19.
429  // These are copied directly from the output of a Python script.
430  Point dx_deta[14] =
431  {
432  -x5/2 + x7/2,
433  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + 3*x5/2 - 3*x7/2 - x9,
434  -x0/2 - x1/2 + 2*x10 - 2*x11 - 2*x12 + x2/2 + x3/2 - 3*x5/2 + 3*x7/2 + 2*x9,
435  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
436  x0/2 + x1/2 + x2/2 + x3/2 - x6 - x8,
437  -x0 - x1 - x2 - x3 + 2*x6 + 2*x8,
438  x0/2 + x1/2 + x2/2 + x3/2 - x6 - x8,
439  x0/4 - x1/4 + x2/4 - x3/4,
440  -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9,
441  x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9,
442  -x0/2 + x1/2 + x2/2 - x3/2 - x6 + x8,
443  x0/2 - x1/2 - x2/2 + x3/2 + x6 - x8,
444  -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
445  x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2
446  };
447 
448  // dx/dxi and dx/deta have 14 components while dx/dzeta has 19.
449  // These are copied directly from the output of a Python script.
450  Point dx_dzeta[19] =
451  {
452  x0/4 + x1/4 + x10 + x11 + x12 + x2/4 + x3/4 - x4 - x5 - x6 - x7 - x8 + x9,
453  -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,
454  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,
455  -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,
456  2*x10 + 2*x11 + 2*x12 - 4*x4 - x5 - x6 - x7 - x8 + 2*x9,
457  x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
458  -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,
459  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,
460  -x0/4 - x1/4 + x10 - x11 - x12 + x2/4 + x3/4 - x5/2 + x7/2 + x9,
461  x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9,
462  -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,
463  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,
464  -x0/4 + x1/4 - x10 - x11 + x12 + x2/4 - x3/4 + x6/2 - x8/2 + x9,
465  -x0/4 + x1/4 - x10 + x11 - x12 - x2/4 + x3/4 + x9,
466  x0/4 - x1/4 + x10 - x11 + x12 + x2/4 - x3/4 - x9,
467  -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2,
468  x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2,
469  -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
470  x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2
471  };
472 
473  // The (xi, eta, zeta) exponents for each of the dx_dxi terms
474  static const int dx_dxi_exponents[14][3] =
475  {
476  {0, 0, 0},
477  {0, 0, 1},
478  {0, 0, 2},
479  {0, 0, 3},
480  {0, 1, 0},
481  {0, 1, 1},
482  {0, 1, 2},
483  {0, 2, 0},
484  {0, 2, 1},
485  {1, 0, 0},
486  {1, 0, 1},
487  {1, 0, 2},
488  {1, 1, 0},
489  {1, 1, 1}
490  };
491 
492  // The (xi, eta, zeta) exponents for each of the dx_deta terms
493  static const int dx_deta_exponents[14][3] =
494  {
495  {0, 0, 0},
496  {0, 0, 1},
497  {0, 0, 2},
498  {0, 0, 3},
499  {0, 1, 0},
500  {0, 1, 1},
501  {0, 1, 2},
502  {1, 0, 0},
503  {1, 0, 1},
504  {1, 0, 2},
505  {1, 1, 0},
506  {1, 1, 1},
507  {2, 0, 0},
508  {2, 0, 1}
509  };
510 
511  // The (xi, eta, zeta) exponents for each of the dx_dzeta terms
512  static const int dx_dzeta_exponents[19][3] =
513  {
514  {0, 0, 0},
515  {0, 0, 1},
516  {0, 0, 2},
517  {0, 0, 3},
518  {0, 0, 4},
519  {0, 1, 0},
520  {0, 1, 1},
521  {0, 1, 2},
522  {0, 1, 3},
523  {1, 0, 0},
524  {1, 0, 1},
525  {1, 0, 2},
526  {1, 0, 3},
527  {1, 1, 0},
528  {1, 1, 1},
529  {1, 2, 0},
530  {1, 2, 1},
531  {2, 1, 0},
532  {2, 1, 1},
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 = Real(-7.1805574131988893873307823958101e-01L),
542  a2 = Real(-5.0580870785392503961340276902425e-01L),
543  a3 = Real(-2.2850430565396735359574351631314e-01L),
544  // Parameters used for zeta quadrature points.
545  b1 = Real( 7.2994024073149732155837979012003e-02L),
546  b2 = Real( 3.4700376603835188472176354340395e-01L),
547  b3 = Real( 7.0500220988849838312239847758405e-01L),
548  // There are 9 unique weight values since there are three
549  // for each of the three unique zeta values.
550  w1 = Real( 4.8498876871878584357513834016440e-02L),
551  w2 = Real( 4.5137737425884574692441981593901e-02L),
552  w3 = Real( 9.2440441384508327195915094925393e-03L),
553  w4 = Real( 7.7598202995005734972022134426305e-02L),
554  w5 = Real( 7.2220379881415319507907170550242e-02L),
555  w6 = Real( 1.4790470621521332351346415188063e-02L),
556  w7 = Real( 1.2415712479200917595523541508209e-01L),
557  w8 = Real( 1.1555260781026451121265147288039e-01L),
558  w9 = Real( 2.3664752994434131762154264300901e-02L);
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<14; ++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<19; ++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 Pyramid13::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 Pyramid13::flip(BoundaryInfo * boundary_info)
721 {
722  swap2nodes(0,1);
723  swap2nodes(2,3);
724  swap2nodes(6,8);
725  swap2nodes(9,10);
726  swap2nodes(11,12);
727  swap2neighbors(1,3);
728  swap2boundarysides(1,3,boundary_info);
729  swap2boundaryedges(1,3,boundary_info);
730  swap2boundaryedges(4,5,boundary_info);
731  swap2boundaryedges(6,7,boundary_info);
732 }
733 
734 
735 ElemType Pyramid13::side_type (const unsigned int s) const
736 {
737  libmesh_assert_less (s, 5);
738  if (s < 4)
739  return TRI6;
740  return QUAD8;
741 }
742 
743 
744 } // 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:3171
static const unsigned int edge_nodes_map[num_edges][nodes_per_edge]
This maps the node of the edge to element node numbers.
static const unsigned int side_nodes_map[num_sides][nodes_per_side]
This maps the node of the side to element node numbers.
virtual std::vector< unsigned int > nodes_on_side(const unsigned int s) const override
virtual bool has_affine_map() 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:2087
virtual unsigned int n_nodes() const override
virtual Real volume() const override
Specialization for computing the volume of a Pyramid13.
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
virtual unsigned int n_second_order_adjacent_vertices(const unsigned int n) const override
IOPackage
libMesh interfaces with several different software packages for the purposes of creating, reading, and writing mesh files.
static const int nodes_per_edge
virtual std::unique_ptr< Elem > build_edge_ptr(const unsigned int i) override
Builds a EDGE3 coincident with edge i.
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:3155
unsigned int p_level() const
Definition: elem.h:2945
virtual unsigned int n_sides() const override
Definition: cell_pyramid.h:83
static const int nodes_per_side
The libMesh namespace provides an interface to certain functionality in the library.
static const int num_children
virtual unsigned int n_vertices() const override
Definition: cell_pyramid.h:88
virtual bool is_edge(const unsigned int i) const override
virtual bool is_vertex(const unsigned int i) const override
T triple_product(const TypeVector< T > &a, const TypeVector< T > &b, const TypeVector< T > &c)
Definition: type_vector.h:1068
void swap4nodes(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four node_ptrs, "rotating" them.
Definition: elem.h:1984
ElemMappingType mapping_type() const
Definition: elem.h:2957
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs.
Definition: elem.h:1933
virtual unsigned int n_sub_elem() const override
FIXME: we don&#39;t yet have a refinement pattern for pyramids...
static const int num_sides
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
Definition: boundary_info.h:57
libmesh_assert(ctx)
static const int num_nodes
Geometric constants for Pyramid13.
virtual void flip(BoundaryInfo *) override final
Flips the element (by swapping node and neighbor pointers) to have a mapping Jacobian of opposite sig...
void swap2neighbors(unsigned int n1, unsigned int n2)
Swaps two neighbor_ptrs.
Definition: elem.h:1943
virtual unsigned int local_side_node(unsigned int side, unsigned int side_node) const override
virtual bool is_face(const unsigned int i) const override
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 override
subdomain_id_type subdomain_id() const
Definition: elem.h:2391
virtual unsigned int n_edges() const override
Definition: cell_pyramid.h:93
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 num_edges
const Node * node_ptr(const unsigned int i) const
Definition: elem.h:2331
ElemType side_type(const unsigned int s) const override final
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=false) override
Builds a QUAD8 or TRI6 coincident with face i.
virtual Real volume() const
Definition: elem.C:3050
void swap4neighbors(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four neighbor_ptrs, "rotating" them.
Definition: elem.h:1994
virtual Order default_order() const override
virtual unsigned int local_edge_node(unsigned int edge, unsigned int edge_node) const override
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
const Point & point(const unsigned int i) const
Definition: elem.h:2277
virtual bool is_node_on_edge(const unsigned int n, const unsigned int e) const override