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"    41     {0, 1, 4, 5, 10,  9, 99, 99, 99}, 
    42     {1, 2, 4, 6, 11, 10, 99, 99, 99}, 
    43     {2, 3, 4, 7, 12, 11, 99, 99, 99}, 
    44     {3, 0, 4, 8,  9, 12, 99, 99, 99}, 
    45     {0, 3, 2, 1,  8,  7,  6,  5, 13}  
    93                                 const unsigned int s)
 const    95   libmesh_assert_less (s, 
n_sides());
   101 std::vector<unsigned>
   104   libmesh_assert_less(s, 
n_sides());
   105   auto trim = (s == 4) ? 0 : 3;
   109 std::vector<unsigned>
   112   libmesh_assert_less(e, 
n_edges());
   117                                 const unsigned int e)
 const   119   libmesh_assert_less (e, 
n_edges());
   143   libmesh_assert_less (s, this->
n_sides());
   157       libmesh_error_msg(
"Invalid side s = " << s);
   164                                         unsigned int side_node)
 const   166   libmesh_assert_less (side, this->
n_sides());
   180                                         unsigned int edge_node)
 const   182   libmesh_assert_less(edge, this->
n_edges());
   192   libmesh_assert_less (i, this->
n_sides());
   194   std::unique_ptr<Elem> face;
   203         face = std::make_unique<Tri6>();
   208         face = std::make_unique<Quad9>();
   212       libmesh_error_msg(
"Invalid side i = " << i);
   216   for (
auto n : face->node_index_range())
   219   face->set_interior_parent(
this);
   220   face->inherit_data_from(*
this);
   228                                 const unsigned int i)
   230   libmesh_assert_less (i, this->
n_sides());
   239         if (!side.get() || side->type() != 
TRI6)
   248         if (!side.get() || side->type() != 
QUAD9)
   256       libmesh_error_msg(
"Invalid side i = " << i);
   259   side->inherit_data_from(*
this);
   262   for (
auto n : side->node_index_range())
   270   return this->simple_build_edge_ptr<Edge3,Pyramid14>(i);
   277   this->simple_build_edge_ptr<Pyramid14>(edge, i, 
EDGE3);
   284                              std::vector<dof_id_type> & )
 const   295         libmesh_not_implemented();
   301         libmesh_not_implemented();
   305       libmesh_error_msg(
"Unsupported IO package " << iop);
   329       libmesh_error_msg(
"Invalid node n = " << n);
   335                                                             const unsigned int v)
 const   337   libmesh_assert_greater_equal (n, this->
n_vertices());
   338   libmesh_assert_less (n, this->
n_nodes());
   351         libmesh_assert_less (v, 2);
   357         unsigned short node_list[8][2] =
   369         return node_list[n-5][v];
   375         libmesh_assert_less (v, 4);
   379         return cast_int<unsigned short>(v);
   383       libmesh_error_msg(
"Invalid n = " << n);
   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,
   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
   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,
   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
   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
   471   static const int dx_dxi_exponents[15][3] =
   491   static const int dx_deta_exponents[15][3] =
   511   static const int dx_dzeta_exponents[20][3] =
   541     a1 = -7.1805574131988893873307823958101e-01_R,
   542     a2 = -5.0580870785392503961340276902425e-01_R,
   543     a3 = -2.2850430565396735359574351631314e-01_R,
   545     b1 = 7.2994024073149732155837979012003e-02_R,
   546     b2 = 3.4700376603835188472176354340395e-01_R,
   547     b3 = 7.0500220988849838312239847758405e-01_R,
   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;
   561   static const Real xi[N][3] =
   592   static const Real eta[N][3] =
   623   static const Real zeta[N][5] =
   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}
   654   static const Real w[N] = {w1, w2, w3, w4, w5, w6, 
   655                             w1, w2, w3, w4, w5, w6, 
   656                             w7, w8, w9, w4, w5, w6, 
   657                             w1, w2, w3, w4, w5, w6, 
   661   for (
int q=0; q<N; ++q)
   665         den2 = (1. - zeta[q][1])*(1. - zeta[q][1]),
   666         den3 = den2*(1. - zeta[q][1]);
   669       Point dx_dxi_q, dx_deta_q;
   670       for (
int c=0; c<15; ++c)
   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];
   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];
   685       for (
int c=0; c<20; ++c)
   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];
   708   libmesh_assert_less (perm_num, 4);
   710   for (
unsigned int i = 0; i != perm_num; ++i)
   746   libmesh_assert_less (s, 5);
 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. 
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. 
Node ** _nodes
Pointers to the nodes we are connected to. 
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
virtual bool is_node_on_side(const unsigned int n, const unsigned int s) const override
static const int num_edges
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
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. 
virtual unsigned int n_sides() const override
The libMesh namespace provides an interface to certain functionality in the library. 
virtual unsigned int n_vertices() const override
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'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)
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. 
ElemMappingType mapping_type() const
void swap2nodes(unsigned int n1, unsigned int n2)
Swaps two node_ptrs. 
The BoundaryInfo class contains information relevant to boundary conditions including storing faces...
static const int num_sides
Geometric constants for all Pyramids. 
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. 
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
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
void swap4neighbors(unsigned int n1, unsigned int n2, unsigned int n3, unsigned int n4)
Swaps four neighbor_ptrs, "rotating" them. 
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)
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. 
dof_id_type node_id(const unsigned int i) const
const Point & point(const unsigned int i) const
static const int nodes_per_edge
virtual std::vector< unsigned int > nodes_on_edge(const unsigned int e) const override
virtual void connectivity(const unsigned int sc, const IOPackage iop, std::vector< dof_id_type > &conn) const override