20 #include "libmesh/fe.h"    21 #include "libmesh/elem.h"    22 #include "libmesh/enum_to_string.h"    32                                       const bool add_p_level)
    37   const Order totalorder = order + add_p_level*elem->
p_level();
    38   libmesh_assert_less(i, n_dofs(elem->
type(), totalorder));
    43   const Real eta  = p(1);
    44   const Real zeta = p(2);
    59               libmesh_assert_less_equal ( std::fabs(xi),   1.0+10*
TOLERANCE );
    60               libmesh_assert_less_equal ( std::fabs(eta),  1.0+10*
TOLERANCE );
    61               libmesh_assert_less_equal ( std::fabs(zeta), 1.0+10*
TOLERANCE );
    66                   return sign * 
RealGradient( -0.125*(1.0-eta-zeta+eta*zeta), 0.0, 0.0 );
    68                   return sign * 
RealGradient( 0.0, -0.125*(1.0+xi-zeta-xi*zeta), 0.0 );
    70                   return sign * 
RealGradient( 0.125*(1.0+eta-zeta-eta*zeta), 0.0, 0.0 );
    72                   return sign * 
RealGradient( 0.0, -0.125*(1.0-xi-zeta+xi*zeta), 0.0 );
    74                   return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0-xi-eta+xi*eta) );
    76                   return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0+xi-eta-xi*eta) );
    78                   return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0+xi+eta+xi*eta) );
    80                   return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0-xi+eta-xi*eta) );
    82                   return sign * 
RealGradient( -0.125*(1.0-eta+zeta-eta*zeta), 0.0, 0.0 );
    84                   return sign * 
RealGradient( 0.0, -0.125*(1.0+xi+zeta+xi*zeta), 0.0 );
    86                   return sign * 
RealGradient( 0.125*(1.0+eta+zeta+eta*zeta), 0.0, 0.0 );
    88                   return sign * 
RealGradient( 0.0, -0.125*(1.0-xi+zeta-xi*zeta), 0.0 );
    91                   libmesh_error_msg(
"Invalid i = " << i);
   107                   return sign * 
RealGradient( -zeta, -zeta, -1.0+xi+eta );
   114                   libmesh_error_msg(
"Invalid i = " << i);
   125       libmesh_error_msg(
"ERROR: Unsupported 3D FE order!: " << totalorder);
   128 #else // LIBMESH_DIM != 3   130   libmesh_not_implemented();
   141   libmesh_error_msg(
"Nedelec elements require the element type \nbecause edge orientation is needed.");
   149                                       const unsigned int i,
   151                                       const bool add_p_level)
   160                                             const unsigned int i,
   161                                             const unsigned int j,
   163                                             const bool add_p_level)
   167   libmesh_assert_less (j, 3);
   169   const Order totalorder = order + add_p_level*elem->
p_level();
   170   libmesh_assert_less(i, n_dofs(elem->
type(), totalorder));
   174   const Real xi   = p(0);
   175   const Real eta  = p(1);
   176   const Real zeta = p(2);
   183         switch (elem->
type())
   191               libmesh_assert_less_equal ( std::fabs(xi),   1.0+10*
TOLERANCE );
   192               libmesh_assert_less_equal ( std::fabs(eta),  1.0+10*
TOLERANCE );
   193               libmesh_assert_less_equal ( std::fabs(zeta), 1.0+10*
TOLERANCE );
   212                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(-1.0+eta) );
   214                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0-eta) );
   216                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0+eta) );
   218                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(-1.0-eta) );
   220                         return sign * 
RealGradient( 0.0, -0.125*(1.0+zeta), 0.0 );
   222                         return sign * 
RealGradient( 0.0, -0.125*(-1.0-zeta), 0.0 );
   225                         libmesh_error_msg(
"Invalid i = " << i);
   241                         return sign * 
RealGradient( -0.125*(-1.0+zeta), 0.0, 0.0 );
   243                         return sign * 
RealGradient( 0.125*(1.0-zeta), 0.0, 0.0 );
   245                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(-1.0+xi) );
   247                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(-1.0-xi) );
   249                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0+xi) );
   251                         return sign * 
RealGradient( 0.0, 0.0, -0.125*(1.0-xi) );
   253                         return sign * 
RealGradient( -0.125*(-1.0-zeta), 0.0, 0.0 );
   255                         return sign * 
RealGradient( 0.125*(1.0+zeta), 0.0, 0.0 );
   258                         libmesh_error_msg(
"Invalid i = " << i);
   274                         return sign * 
RealGradient( -0.125*(-1.0+eta), 0.0, 0.0 );
   276                         return sign * 
RealGradient( 0.0, -0.125*(-1.0-xi), 0.0 );
   278                         return sign * 
RealGradient( 0.125*(-1.0-eta), 0.0, 0.0 );
   280                         return sign * 
RealGradient( 0.0, -0.125*(-1.0+xi), 0.0 );
   282                         return sign * 
RealGradient( -0.125*(1.0-eta), 0.0, 0.0 );
   284                         return sign * 
RealGradient( 0.0, -0.125*(1.0+xi), 0.0 );
   286                         return sign * 
RealGradient( 0.125*(1.0+eta), 0.0, 0.0 );
   288                         return sign * 
RealGradient( 0.0, -0.125*(1.0-xi), 0.0 );
   291                         libmesh_error_msg(
"Invalid i = " << i);
   297                   libmesh_error_msg(
"Invalid j = " << j);
   325                         libmesh_error_msg(
"Invalid i = " << i);
   348                         libmesh_error_msg(
"Invalid i = " << i);
   371                           libmesh_error_msg(
"Invalid i = " << i);
   377                   libmesh_error_msg(
"Invalid j = " << j);
   387       libmesh_error_msg(
"ERROR: Unsupported 3D FE order!: " << totalorder);
   390 #else // LIBMESH_DIM != 3   392   libmesh_not_implemented();
   404   libmesh_error_msg(
"Nedelec elements require the element type \nbecause edge orientation is needed.");
   412                                             const unsigned int i,
   413                                             const unsigned int j,
   415                                             const bool add_p_level)
   421 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES   426                                                    const unsigned int i,
   427                                                    const unsigned int j,
   429                                                    const bool add_p_level)
   441   libmesh_assert_less (j, 6);
   443   const Order totalorder = order + add_p_level*elem->
p_level();
   444   libmesh_assert_less(i, n_dofs(elem->
type(), totalorder));
   453         switch (elem->
type())
   488                         libmesh_error_msg(
"Invalid i = " << i);
   515                         libmesh_error_msg(
"Invalid i = " << i);
   542                         libmesh_error_msg(
"Invalid i = " << i);
   548                   libmesh_error_msg(
"Invalid j = " << j);
   566       libmesh_error_msg(
"ERROR: Unsupported 3D FE order!: " << totalorder);
   569 #else // LIBMESH_DIM != 3   572   libmesh_not_implemented();
   584   libmesh_error_msg(
"Nedelec elements require the element type \nbecause edge orientation is needed.");
   592                                                    const unsigned int i,
   593                                                    const unsigned int j,
   595                                                    const bool add_p_level)
 class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
ElemType
Defines an enum for geometric element types. 
RealVectorValue RealGradient
Order
defines an enum for polynomial orders. 
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
static constexpr Real TOLERANCE
This is the base class from which all geometric element types are derived. 
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
unsigned int p_level() const
OrderWrapper order
The approximation order of the element. 
The libMesh namespace provides an interface to certain functionality in the library. 
void libmesh_ignore(const Args &...)
bool positive_edge_orientation(const unsigned int i) const
std::string enum_to_string(const T e)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space. 
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)