25 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 36 const Order libmesh_dbg_var(order),
41 const Real xi2 = xi*xi;
46 libmesh_assert_less_equal (order,
SEVENTH);
61 case 0:
return 1./2.-1./2.*xi;
62 case 1:
return 1./2.+1./2.*xi;
63 case 2:
return 1./4. *2.4494897427831780982*(xi2-1.);
64 case 3:
return 1./4. *3.1622776601683793320*(xi2-1.)*xi;
65 case 4:
return 1./16. *3.7416573867739413856*((5.*xi2-6.)*xi2+1.);
66 case 5:
return 3./16. *1.4142135623730950488*(3.+(-10.+7.*xi2)*xi2)*xi;
67 case 6:
return 1./32. *4.6904157598234295546*(-1.+(15.+(-35.+21.*xi2)*xi2)*xi2);
68 case 7:
return 1./32. *5.0990195135927848300*(-5.+(35.+(-63.+33.*xi2)*xi2)*xi2)*xi;
69 case 8:
return 1./256.*5.4772255750516611346*(5.+(-140.+(630.+(-924.+429.*xi2)*xi2)*xi2)*xi2);
72 libmesh_error_msg(
"Invalid shape function index!");
93 const Order libmesh_dbg_var(order),
95 const unsigned int libmesh_dbg_var(j),
99 libmesh_assert_equal_to (j, 0);
101 const Real xi = p(0);
102 const Real xi2 = xi*xi;
106 libmesh_assert_less_equal (order,
SEVENTH);
120 case 0:
return -1./2.;
122 case 2:
return 1./2.*2.4494897427831780982*xi;
123 case 3:
return -1./4.*3.1622776601683793320+3./4.*3.1622776601683793320*xi2;
124 case 4:
return 1./16.*3.7416573867739413856*(-12.+20*xi2)*xi;
125 case 5:
return 9./16.*1.4142135623730950488+(-45./8.*1.4142135623730950488+105./16.*1.4142135623730950488*xi2)*xi2;
126 case 6:
return 1./32.*4.6904157598234295546*(30.+(-140.+126.*xi2)*xi2)*xi;
127 case 7:
return -5./32.*5.0990195135927848300+(105./32.*5.0990195135927848300+(-315./32.*5.0990195135927848300+231./32.*5.0990195135927848300*xi2)*xi2)*xi2;
128 case 8:
return 1./256.*5.4772255750516611346*(-280.+(2520.+(-5544.+3432.*xi2)*xi2)*xi2)*xi;
131 libmesh_error_msg(
"Invalid shape function index!");
140 const unsigned int i,
141 const unsigned int j,
144 libmesh_assert(elem);
159 static bool warning_given =
false;
162 libMesh::err <<
"Second derivatives for Szabab elements " 163 <<
" are not yet implemented!" 166 warning_given =
true;
179 static bool warning_given =
false;
182 libMesh::err <<
"Second derivatives for Szabab elements " 183 <<
" are not yet implemented!" 186 warning_given =
true;
193 #endif //LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
The base class for all geometric element types.
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
OStreamProxy err(std::cerr)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual ElemType type() const =0
A geometric point in (x,y,z) space.
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)