20 #include "libmesh/side.h" 
   21 #include "libmesh/cell_pyramid14.h" 
   22 #include "libmesh/edge_edge3.h" 
   23 #include "libmesh/face_tri6.h" 
   24 #include "libmesh/face_quad9.h" 
   25 #include "libmesh/enum_io_package.h" 
   26 #include "libmesh/enum_order.h" 
   45     {0, 1, 4, 5, 10,  9, 99, 99, 99}, 
 
   46     {1, 2, 4, 6, 11, 10, 99, 99, 99}, 
 
   47     {2, 3, 4, 7, 12, 11, 99, 99, 99}, 
 
   48     {3, 0, 4, 8,  9, 12, 99, 99, 99}, 
 
   49     {0, 3, 2, 1,  8,  7,  6,  5, 13}  
 
   99                                 const unsigned int s)
 const 
  101   libmesh_assert_less (s, 
n_sides());
 
  107 std::vector<unsigned>
 
  110   libmesh_assert_less(s, 
n_sides());
 
  111   auto trim = (s == 4) ? 0 : 3;
 
  116                                 const unsigned int e)
 const 
  118   libmesh_assert_less (e, 
n_edges());
 
  144   libmesh_assert_less (s, this->
n_sides());
 
  158       libmesh_error_msg(
"Invalid side s = " << s);
 
  165                                         unsigned int side_node)
 const 
  167   libmesh_assert_less (side, this->
n_sides());
 
  170   libmesh_assert_less(side_node, 9);
 
  182   libmesh_assert_less (i, this->
n_sides());
 
  192           return libmesh_make_unique<Side<Tri6,Pyramid14>>(
this,i);
 
  195           return libmesh_make_unique<Side<Quad9,Pyramid14>>(
this,i);
 
  198           libmesh_error_msg(
"Invalid side i = " << i);
 
  205       std::unique_ptr<Elem> face;
 
  214             face = libmesh_make_unique<Tri6>();
 
  219             face = libmesh_make_unique<Quad9>();
 
  223           libmesh_error_msg(
"Invalid side i = " << i);
 
  229       for (
auto n : face->node_index_range())
 
  230         face->set_node(n) = this->
node_ptr(Pyramid14::side_nodes_map[i][n]);
 
  239                                 const unsigned int i)
 
  241   libmesh_assert_less (i, this->
n_sides());
 
  250         if (!side.get() || side->type() != 
TRI6)
 
  259         if (!side.get() || side->type() != 
QUAD9)
 
  267       libmesh_error_msg(
"Invalid side i = " << i);
 
  273   for (
auto n : side->node_index_range())
 
  274     side->set_node(n) = this->
node_ptr(Pyramid14::side_nodes_map[i][n]);
 
  281   libmesh_assert_less (i, this->
n_edges());
 
  283   return libmesh_make_unique<SideEdge<Edge3,Pyramid14>>(
this,i);
 
  290                              std::vector<dof_id_type> & )
 const 
  301         libmesh_not_implemented();
 
  307         libmesh_not_implemented();
 
  311       libmesh_error_msg(
"Unsupported IO package " << iop);
 
  335       libmesh_error_msg(
"Invalid node n = " << n);
 
  341                                                             const unsigned int v)
 const 
  343   libmesh_assert_greater_equal (n, this->
n_vertices());
 
  344   libmesh_assert_less (n, this->
n_nodes());
 
  357         libmesh_assert_less (v, 2);
 
  363         unsigned short node_list[8][2] =
 
  375         return node_list[n-5][v];
 
  381         libmesh_assert_less (v, 4);
 
  385         return cast_int<unsigned short>(v);
 
  389       libmesh_error_msg(
"Invalid n = " << n);
 
  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,
 
  424       -x0/2 - x1/2 + x2/2 + x3/2 + x5 - x7,
 
  425       x0/2 + x1/2 - x2/2 - x3/2 - x5 + x7,
 
  426       x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
 
  434       x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + 3*x5/2 - 3*x7/2 - x9,
 
  435       -x0/2 - x1/2 + 2*x10 - 2*x11 - 2*x12 + x2/2 + x3/2 - 3*x5/2 + 3*x7/2 + 2*x9,
 
  436       x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
 
  440       x0/4 - x1/4 + x2/4 - x3/4,
 
  441       -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9,
 
  442       x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9,
 
  443       -x0/2 + x1/2 + x2/2 - x3/2 - x6 + x8,
 
  444       x0/2 - x1/2 - x2/2 + x3/2 + x6 - x8,
 
  445       -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
 
  446       x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2,
 
  447       x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
 
  454       -x0/4 - x1/4 + x10 + x11 + x12 - 2*x13 - x2/4 - x3/4 - x4 + x9,
 
  455       5*x0/4 + 5*x1/4 - 5*x10 - 5*x11 - 5*x12 + 8*x13 + 5*x2/4 + 5*x3/4 + 7*x4 - 5*x9,
 
  456       -9*x0/4 - 9*x1/4 + 9*x10 + 9*x11 + 9*x12 - 12*x13 - 9*x2/4 - 9*x3/4 - 15*x4 + 9*x9,
 
  457       7*x0/4 + 7*x1/4 - 7*x10 - 7*x11 - 7*x12 + 8*x13 + 7*x2/4 + 7*x3/4 + 13*x4 - 7*x9,
 
  458       -x0/2 - x1/2 + 2*x10 + 2*x11 + 2*x12 - 2*x13 - x2/2 - x3/2 - 4*x4 + 2*x9,
 
  459       x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9,
 
  460       -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,
 
  461       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,
 
  462       -x0/4 - x1/4 + x10 - x11 - x12 + x2/4 + x3/4 - x5/2 + x7/2 + x9,
 
  463       x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9,
 
  464       -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,
 
  465       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,
 
  466       -x0/4 + x1/4 - x10 - x11 + x12 + x2/4 - x3/4 + x6/2 - x8/2 + x9,
 
  467       -x0/4 + x1/4 - x10 + x11 - x12 - x2/4 + x3/4 + x9,
 
  468       x0/4 - x1/4 + x10 - x11 + x12 + x2/4 - x3/4 - x9,
 
  469       -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2,
 
  470       x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2,
 
  471       -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2,
 
  472       x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2,
 
  473       x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8
 
  477   static const int dx_dxi_exponents[15][3] =
 
  497   static const int dx_deta_exponents[15][3] =
 
  517   static const int dx_dzeta_exponents[20][3] =
 
  547     a1 = 
Real(-7.1805574131988893873307823958101e-01L),
 
  548     a2 = 
Real(-5.0580870785392503961340276902425e-01L),
 
  549     a3 = 
Real(-2.2850430565396735359574351631314e-01L),
 
  551     b1 = 
Real( 7.2994024073149732155837979012003e-02L),
 
  552     b2 = 
Real( 3.4700376603835188472176354340395e-01L),
 
  553     b3 = 
Real( 7.0500220988849838312239847758405e-01L),
 
  556     w1 = 
Real( 4.8498876871878584357513834016440e-02L),
 
  557     w2 = 
Real( 4.5137737425884574692441981593901e-02L),
 
  558     w3 = 
Real( 9.2440441384508327195915094925393e-03L),
 
  559     w4 = 
Real( 7.7598202995005734972022134426305e-02L),
 
  560     w5 = 
Real( 7.2220379881415319507907170550242e-02L),
 
  561     w6 = 
Real( 1.4790470621521332351346415188063e-02L),
 
  562     w7 = 
Real( 1.2415712479200917595523541508209e-01L),
 
  563     w8 = 
Real( 1.1555260781026451121265147288039e-01L),
 
  564     w9 = 
Real( 2.3664752994434131762154264300901e-02L);
 
  567   static const Real xi[N][3] =
 
  598   static const Real eta[N][3] =
 
  629   static const Real zeta[N][5] =
 
  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       { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
 
  653       { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
 
  654       { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3},
 
  655       { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1},
 
  656       { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2},
 
  657       { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}
 
  660   static const Real w[N] = {w1, w2, w3, w4, w5, w6, 
 
  661                             w1, w2, w3, w4, w5, w6, 
 
  662                             w7, w8, w9, w4, w5, w6, 
 
  663                             w1, w2, w3, w4, w5, w6, 
 
  667   for (
int q=0; q<N; ++q)
 
  671         den2 = (1. - zeta[q][1])*(1. - zeta[q][1]),
 
  672         den3 = den2*(1. - zeta[q][1]);
 
  675       Point dx_dxi_q, dx_deta_q;
 
  676       for (
int c=0; c<15; ++c)
 
  679             xi[q][dx_dxi_exponents[c][0]]*
 
  680             eta[q][dx_dxi_exponents[c][1]]*
 
  681             zeta[q][dx_dxi_exponents[c][2]]*dx_dxi[c];
 
  684             xi[q][dx_deta_exponents[c][0]]*
 
  685             eta[q][dx_deta_exponents[c][1]]*
 
  686             zeta[q][dx_deta_exponents[c][2]]*dx_deta[c];
 
  691       for (
int c=0; c<20; ++c)
 
  694             xi[q][dx_dzeta_exponents[c][0]]*
 
  695             eta[q][dx_dzeta_exponents[c][1]]*
 
  696             zeta[q][dx_dzeta_exponents[c][2]]*dx_dzeta[c];