20 #include "libmesh/fe.h" 21 #include "libmesh/elem.h" 22 #include "libmesh/fe_lagrange_shape_1D.h" 23 #include "libmesh/enum_to_string.h" 24 #include "libmesh/cell_c0polyhedron.h" 25 #include "libmesh/tensor_value.h" 49 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 52 Real fe_lagrange_3D_shape_second_deriv(
const ElemType type,
59 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 76 const std::vector<Point> & p,
77 std::vector<std::vector<OutputShape>> & v,
78 const bool add_p_level)
86 (elem,o,p,v,add_p_level);
92 const unsigned int n_sf = v.size();
106 libmesh_assert_less_equal (n_sf, 8);
109 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0};
110 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1};
111 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1};
115 const Point & q_point = p[qp];
117 const Real xi = q_point(0);
118 const Real eta = q_point(1);
119 const Real zeta = q_point(2);
122 Real one_d_shapes[3][2] = {
131 v[i][qp] = one_d_shapes[0][i0[i]] *
132 one_d_shapes[1][i1[i]] *
133 one_d_shapes[2][i2[i]];
154 libmesh_fallthrough();
157 libmesh_assert_less_equal (n_sf, 27);
163 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 0, 2, 2, 1, 2, 0, 2, 2};
164 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 2, 0, 2, 1, 2, 2, 2};
165 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1, 0, 0, 0, 0, 2, 2, 2, 2, 1, 1, 1, 1, 0, 2, 2, 2, 2, 1, 2};
169 const Point & q_point = p[qp];
171 const Real xi = q_point(0);
172 const Real eta = q_point(1);
173 const Real zeta = q_point(2);
176 Real one_d_shapes[3][3] = {
188 v[i][qp] = one_d_shapes[0][i0[i]] *
189 one_d_shapes[1][i1[i]] *
190 one_d_shapes[2][i2[i]];
202 libmesh_error_msg(
"ERROR: Unsupported 3D FE order on HEX!: " << o);
204 #else // LIBMESH_DIM != 3 206 libmesh_not_implemented();
207 #endif // LIBMESH_DIM == 3 214 const unsigned int i,
215 const std::vector<Point> & p,
216 std::vector<OutputShape> & v,
217 const bool add_p_level)
220 (elem,o,i,p,v,add_p_level);
227 const unsigned int i,
228 const unsigned int j,
229 const std::vector<Point> & p,
230 std::vector<OutputShape> & v,
231 const bool add_p_level)
234 (elem,o,i,j,p,v,add_p_level);
241 const std::vector<Point> & p,
242 std::vector<std::vector<OutputShape>> * comps[3],
243 const bool add_p_level)
251 (elem,o,p,comps,add_p_level);
260 const unsigned int n_sf = comps[0]->size();
274 libmesh_assert_equal_to (n_sf, 8);
277 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0};
278 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1};
279 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1};
283 const Point & q_point = p[qp];
285 const Real xi = q_point(0);
286 const Real eta = q_point(1);
287 const Real zeta = q_point(2);
290 Real one_d_shapes[3][2] = {
299 Real one_d_derivs[3][2] = {
309 (*comps[0])[i][qp] = one_d_derivs[0][i0[i]] *
310 one_d_shapes[1][i1[i]] *
311 one_d_shapes[2][i2[i]];
312 (*comps[1])[i][qp] = one_d_shapes[0][i0[i]] *
313 one_d_derivs[1][i1[i]] *
314 one_d_shapes[2][i2[i]];
315 (*comps[2])[i][qp] = one_d_shapes[0][i0[i]] *
316 one_d_shapes[1][i1[i]] *
317 one_d_derivs[2][i2[i]];
339 libmesh_fallthrough();
342 libmesh_assert_less_equal (n_sf, 27);
348 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 0, 2, 2, 1, 2, 0, 2, 2};
349 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 2, 0, 2, 1, 2, 2, 2};
350 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1, 0, 0, 0, 0, 2, 2, 2, 2, 1, 1, 1, 1, 0, 2, 2, 2, 2, 1, 2};
354 const Point & q_point = p[qp];
356 const Real xi = q_point(0);
357 const Real eta = q_point(1);
358 const Real zeta = q_point(2);
361 Real one_d_shapes[3][3] = {
373 Real one_d_derivs[3][3] = {
386 (*comps[0])[i][qp] = one_d_derivs[0][i0[i]] *
387 one_d_shapes[1][i1[i]] *
388 one_d_shapes[2][i2[i]];
389 (*comps[1])[i][qp] = one_d_shapes[0][i0[i]] *
390 one_d_derivs[1][i1[i]] *
391 one_d_shapes[2][i2[i]];
392 (*comps[2])[i][qp] = one_d_shapes[0][i0[i]] *
393 one_d_shapes[1][i1[i]] *
394 one_d_derivs[2][i2[i]];
407 libmesh_error_msg(
"ERROR: Unsupported 3D FE order on HEX!: " << o);
409 #else // LIBMESH_DIM != 3 411 libmesh_not_implemented();
412 #endif // LIBMESH_DIM == 3 421 const unsigned int i,
424 return fe_lagrange_3D_shape<LAGRANGE>(type, order,
nullptr, i, p);
432 const unsigned int i,
435 return fe_lagrange_3D_shape<L2_LAGRANGE>(type, order,
nullptr, i, p);
443 const unsigned int i,
445 const bool add_p_level)
450 return fe_lagrange_3D_shape<LAGRANGE>(elem->
type(), order + add_p_level*elem->
p_level(), elem, i, p);
458 const unsigned int i,
460 const bool add_p_level)
465 return fe_lagrange_3D_shape<L2_LAGRANGE>(elem->
type(), order + add_p_level*elem->
p_level(), elem, i, p);
473 const unsigned int i,
475 const bool add_p_level)
478 return fe_lagrange_3D_shape<LAGRANGE>(elem->
type(), fet.
order + add_p_level*elem->
p_level(), elem, i, p);
486 const unsigned int i,
488 const bool add_p_level)
491 return fe_lagrange_3D_shape<L2_LAGRANGE>(elem->
type(), fet.
order + add_p_level*elem->
p_level(), elem, i, p);
497 const unsigned int i,
498 const unsigned int j,
501 return fe_lagrange_3D_shape_deriv<LAGRANGE>(type, order,
nullptr, i, j, p);
509 const unsigned int i,
510 const unsigned int j,
513 return fe_lagrange_3D_shape_deriv<L2_LAGRANGE>(type, order,
nullptr, i, j, p);
521 const unsigned int i,
522 const unsigned int j,
524 const bool add_p_level)
529 return fe_lagrange_3D_shape_deriv<LAGRANGE>(elem->
type(), order + add_p_level*elem->
p_level(), elem, i, j, p);
536 const unsigned int i,
537 const unsigned int j,
539 const bool add_p_level)
544 return fe_lagrange_3D_shape_deriv<L2_LAGRANGE>(elem->
type(), order + add_p_level*elem->
p_level(), elem, i, j, p);
551 const unsigned int i,
552 const unsigned int j,
554 const bool add_p_level)
557 return fe_lagrange_3D_shape_deriv<LAGRANGE>(elem->
type(), fet.
order + add_p_level*elem->
p_level(), elem, i, j, p);
564 const unsigned int i,
565 const unsigned int j,
567 const bool add_p_level)
570 return fe_lagrange_3D_shape_deriv<L2_LAGRANGE>(elem->
type(), fet.
order + add_p_level*elem->
p_level(), elem, i, j, p);
573 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 578 const unsigned int i,
579 const unsigned int j,
582 return fe_lagrange_3D_shape_second_deriv<LAGRANGE>(type, order,
nullptr, i, j, p);
590 const unsigned int i,
591 const unsigned int j,
594 return fe_lagrange_3D_shape_second_deriv<L2_LAGRANGE>(type, order,
nullptr, i, j, p);
602 const unsigned int i,
603 const unsigned int j,
605 const bool add_p_level)
610 return fe_lagrange_3D_shape_second_deriv<LAGRANGE>
611 (elem->
type(), order + add_p_level*elem->
p_level(), elem, i, j, p);
619 const unsigned int i,
620 const unsigned int j,
622 const bool add_p_level)
627 return fe_lagrange_3D_shape_second_deriv<L2_LAGRANGE>
628 (elem->
type(), order + add_p_level*elem->
p_level(), elem, i, j, p);
635 const unsigned int i,
636 const unsigned int j,
638 const bool add_p_level)
641 return fe_lagrange_3D_shape_second_deriv<LAGRANGE>
650 const unsigned int i,
651 const unsigned int j,
653 const bool add_p_level)
656 return fe_lagrange_3D_shape_second_deriv<L2_LAGRANGE>
661 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 671 template <FEFamily T>
675 const unsigned int i,
692 libmesh_assert_less (i, 8);
695 const Real xi = p(0);
696 const Real eta = p(1);
697 const Real zeta = p(2);
700 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0};
701 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1};
702 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1};
714 libmesh_assert_less (i, 4);
717 const Real zeta1 = p(0);
718 const Real zeta2 = p(1);
719 const Real zeta3 = p(2);
720 const Real zeta0 = 1. - zeta1 - zeta2 - zeta3;
737 libmesh_error_msg(
"Invalid i = " << i);
748 libmesh_assert_less (i, 6);
753 Point p2d(p(0),p(1));
757 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1};
758 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2};
770 libmesh_assert_less (i, 5);
772 const Real xi = p(0);
773 const Real eta = p(1);
774 const Real zeta = p(2);
775 const Real eps = 1.e-35;
780 return .25*(zeta + xi - 1.)*(zeta + eta - 1.)/((1. - zeta) + eps);
783 return .25*(zeta - xi - 1.)*(zeta + eta - 1.)/((1. - zeta) + eps);
786 return .25*(zeta - xi - 1.)*(zeta - eta - 1.)/((1. - zeta) + eps);
789 return .25*(zeta + xi - 1.)*(zeta - eta - 1.)/((1. - zeta) + eps);
795 libmesh_error_msg(
"Invalid i = " << i);
803 libmesh_error_msg(
"Code (see stack trace) used an outdated FE function overload.\n" 804 "Shape functions on a polyhedron are not defined by ElemType alone.");
808 const C0Polyhedron & poly = *cast_ptr<const C0Polyhedron *>(elem);
816 libmesh_assert_less(s, poly.n_subelements());
818 const auto subtet = poly.subelement(s);
822 if (nodei == subtet[0])
824 if (nodei == subtet[1])
826 if (nodei == subtet[2])
828 if (nodei == subtet[3])
850 libmesh_assert_less (i, 20);
852 const Real xi = p(0);
853 const Real eta = p(1);
854 const Real zeta = p(2);
858 const Real x = .5*(xi + 1.);
859 const Real y = .5*(eta + 1.);
860 const Real z = .5*(zeta + 1.);
865 return (1. - x)*(1. - y)*(1. - z)*(1. - 2.*x - 2.*y - 2.*z);
868 return x*(1. - y)*(1. - z)*(2.*x - 2.*y - 2.*z - 1.);
871 return x*y*(1. - z)*(2.*x + 2.*y - 2.*z - 3.);
874 return (1. - x)*y*(1. - z)*(2.*y - 2.*x - 2.*z - 1.);
877 return (1. - x)*(1. - y)*z*(2.*z - 2.*x - 2.*y - 1.);
880 return x*(1. - y)*z*(2.*x - 2.*y + 2.*z - 3.);
883 return x*y*z*(2.*x + 2.*y + 2.*z - 5.);
886 return (1. - x)*y*z*(2.*y - 2.*x + 2.*z - 3.);
889 return 4.*x*(1. - x)*(1. - y)*(1. - z);
892 return 4.*x*y*(1. - y)*(1. - z);
895 return 4.*x*(1. - x)*y*(1. - z);
898 return 4.*(1. - x)*y*(1. - y)*(1. - z);
901 return 4.*(1. - x)*(1. - y)*z*(1. - z);
904 return 4.*x*(1. - y)*z*(1. - z);
907 return 4.*x*y*z*(1. - z);
910 return 4.*(1. - x)*y*z*(1. - z);
913 return 4.*x*(1. - x)*(1. - y)*z;
916 return 4.*x*y*(1. - y)*z;
919 return 4.*x*(1. - x)*y*z;
922 return 4.*(1. - x)*y*(1. - y)*z;
925 libmesh_error_msg(
"Invalid i = " << i);
932 "High order on first order elements only supported for L2 families");
933 libmesh_fallthrough();
936 libmesh_assert_less (i, 27);
939 const Real xi = p(0);
940 const Real eta = p(1);
941 const Real zeta = p(2);
947 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 0, 2, 2, 1, 2, 0, 2, 2};
948 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 2, 0, 2, 1, 2, 2, 2};
949 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1, 0, 0, 0, 0, 2, 2, 2, 2, 1, 1, 1, 1, 0, 2, 2, 2, 2, 1, 2};
959 "High order on first order elements only supported for L2 families");
960 libmesh_fallthrough();
962 libmesh_assert_less (i, 10);
963 libmesh_fallthrough();
966 libmesh_assert_less (i, 14);
969 const Real zeta1 = p(0);
970 const Real zeta2 = p(1);
971 const Real zeta3 = p(2);
972 const Real zeta0 = 1. - zeta1 - zeta2 - zeta3;
977 return zeta0*(2.*zeta0 - 1.);
980 return zeta1*(2.*zeta1 - 1.);
983 return zeta2*(2.*zeta2 - 1.);
986 return zeta3*(2.*zeta3 - 1.);
989 return 4.*zeta0*zeta1;
992 return 4.*zeta1*zeta2;
995 return 4.*zeta2*zeta0;
998 return 4.*zeta0*zeta3;
1001 return 4.*zeta1*zeta3;
1004 return 4.*zeta2*zeta3;
1007 libmesh_error_msg(
"Invalid i = " << i);
1014 libmesh_assert_less (i, 15);
1016 const Real xi = p(0);
1017 const Real eta = p(1);
1018 const Real zeta = p(2);
1023 return (1. - zeta)*(xi + eta - 1.)*(xi + eta + 0.5*zeta);
1026 return (1. - zeta)*xi*(xi - 1. - 0.5*zeta);
1029 return (1. - zeta)*eta*(eta - 1. - 0.5*zeta);
1032 return (1. + zeta)*(xi + eta - 1.)*(xi + eta - 0.5*zeta);
1035 return (1. + zeta)*xi*(xi - 1. + 0.5*zeta);
1038 return (1. + zeta)*eta*(eta - 1. + 0.5*zeta);
1041 return 2.*(1. - zeta)*xi*(1. - xi - eta);
1044 return 2.*(1. - zeta)*xi*eta;
1047 return 2.*(1. - zeta)*eta*(1. - xi - eta);
1050 return (1. - zeta)*(1. + zeta)*(1. - xi - eta);
1053 return (1. - zeta)*(1. + zeta)*xi;
1056 return (1. - zeta)*(1. + zeta)*eta;
1059 return 2.*(1. + zeta)*xi*(1. - xi - eta);
1062 return 2.*(1. + zeta)*xi*eta;
1065 return 2.*(1. + zeta)*eta*(1. - xi - eta);
1068 libmesh_error_msg(
"Invalid i = " << i);
1075 "High order on first order elements only supported for L2 families");
1076 libmesh_fallthrough();
1081 libmesh_assert_less (i, 18);
1086 Point p2d(p(0),p(1));
1090 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2};
1091 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2, 3, 4, 5, 0, 1, 2, 3, 4, 5, 3, 4, 5};
1102 libmesh_assert_less (i, 13);
1104 const Real xi = p(0);
1105 const Real eta = p(1);
1106 const Real zeta = p(2);
1107 const Real eps = 1.e-35;
1111 Real den = (1. - zeta + eps);
1116 return 0.25*(-xi - eta - 1.)*((1. - xi)*(1. - eta) - zeta + xi*eta*zeta/den);
1119 return 0.25*(-eta + xi - 1.)*((1. + xi)*(1. - eta) - zeta - xi*eta*zeta/den);
1122 return 0.25*(xi + eta - 1.)*((1. + xi)*(1. + eta) - zeta + xi*eta*zeta/den);
1125 return 0.25*(eta - xi - 1.)*((1. - xi)*(1. + eta) - zeta - xi*eta*zeta/den);
1128 return zeta*(2.*zeta - 1.);
1131 return 0.5*(1. + xi - zeta)*(1. - xi - zeta)*(1. - eta - zeta)/den;
1134 return 0.5*(1. + eta - zeta)*(1. - eta - zeta)*(1. + xi - zeta)/den;
1137 return 0.5*(1. + xi - zeta)*(1. - xi - zeta)*(1. + eta - zeta)/den;
1140 return 0.5*(1. + eta - zeta)*(1. - eta - zeta)*(1. - xi - zeta)/den;
1143 return zeta*(1. - xi - zeta)*(1. - eta - zeta)/den;
1146 return zeta*(1. + xi - zeta)*(1. - eta - zeta)/den;
1149 return zeta*(1. + eta - zeta)*(1. + xi - zeta)/den;
1152 return zeta*(1. - xi - zeta)*(1. + eta - zeta)/den;
1155 libmesh_error_msg(
"Invalid i = " << i);
1164 "High order on first order elements only supported for L2 families");
1165 libmesh_fallthrough();
1169 libmesh_assert_less (i, 14);
1171 const Real xi = p(0);
1172 const Real eta = p(1);
1173 const Real zeta = p(2);
1174 const Real eps = 1.e-35;
1179 p1 = 0.5*(1. - eta - zeta),
1180 p2 = 0.5*(1. + xi - zeta),
1181 p3 = 0.5*(1. + eta - zeta),
1182 p4 = 0.5*(1. - xi - zeta);
1187 den = (-1. + zeta + eps),
1193 return p4*p1*(xi*eta - zeta + zeta*zeta)/den2;
1196 return -p1*p2*(xi*eta + zeta - zeta*zeta)/den2;
1199 return p2*p3*(xi*eta - zeta + zeta*zeta)/den2;
1202 return -p3*p4*(xi*eta + zeta - zeta*zeta)/den2;
1205 return zeta*(2.*zeta - 1.);
1208 return -4.*p2*p1*p4*eta/den2;
1211 return 4.*p1*p2*p3*xi/den2;
1214 return 4.*p2*p3*p4*eta/den2;
1217 return -4.*p3*p4*p1*xi/den2;
1220 return -4.*p1*p4*zeta/den;
1223 return -4.*p2*p1*zeta/den;
1226 return -4.*p3*p2*zeta/den;
1229 return -4.*p4*p3*zeta/den;
1232 return 16.*p1*p2*p3*p4/den2;
1235 libmesh_error_msg(
"Invalid i = " << i);
1252 libmesh_assert_less (i, 20);
1260 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2, 0, 1};
1262 Point p2d(p(0),p(1));
1270 const Real bubbleval =
1275 return mainval - bubbleval / 9;
1277 return mainval + bubbleval * (
Real(4) / 9);
1283 libmesh_assert_less (i, 21);
1288 Point p2d(p(0),p(1));
1292 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2, 0, 1, 2};
1293 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2, 3, 4, 5, 0, 1, 2, 3, 4, 5, 3, 4, 5, 6, 6, 6};
1302 libmesh_assert_less (i, 18);
1304 const Real xi = p(0);
1305 const Real eta = p(1);
1306 const Real zeta = p(2);
1307 const Real eps = 1.e-35;
1312 p1 = 0.5*(1. - eta - zeta),
1313 p2 = 0.5*(1. + xi - zeta),
1314 p3 = 0.5*(1. + eta - zeta),
1315 p4 = 0.5*(1. - xi - zeta);
1320 den = (-1. + zeta + eps),
1327 constexpr
Real alpha = 0.5;
1329 bub_f1 = ((1-alpha)*(1-zeta) + alpha*(-eta)),
1330 bub_f2 = ((1-alpha)*(1-zeta) + alpha*(xi)),
1331 bub_f3 = ((1-alpha)*(1-zeta) + alpha*(eta)),
1332 bub_f4 = ((1-alpha)*(1-zeta) + alpha*(-xi));
1335 bub1 = bub_f1*p1*p2*p4*zeta/den2,
1336 bub2 = bub_f2*p1*p2*p3*zeta/den2,
1337 bub3 = bub_f3*p2*p3*p4*zeta/den2,
1338 bub4 = bub_f4*p1*p3*p4*zeta/den2;
1343 return p4*p1*(xi*eta - zeta + zeta*zeta)/den2 + 3*(bub1+bub4);
1346 return -p1*p2*(xi*eta + zeta - zeta*zeta)/den2 + 3*(bub1+bub2);
1349 return p2*p3*(xi*eta - zeta + zeta*zeta)/den2 + 3*(bub2+bub3);
1352 return -p3*p4*(xi*eta + zeta - zeta*zeta)/den2 + 3*(bub3+bub4);
1355 return zeta*(2.*zeta - 1.) + 3*(bub1+bub2+bub3+bub4);
1358 return -4.*p2*p1*p4*eta/den2 - 12*bub1;
1361 return 4.*p1*p2*p3*xi/den2 - 12*bub2;
1364 return 4.*p2*p3*p4*eta/den2 - 12*bub3;
1367 return -4.*p3*p4*p1*xi/den2 - 12*bub4;
1370 return -4.*p1*p4*zeta/den - 12*(bub1+bub4);
1373 return -4.*p2*p1*zeta/den - 12*(bub1+bub2);
1376 return -4.*p3*p2*zeta/den - 12*(bub2+bub3);
1379 return -4.*p4*p3*zeta/den - 12*(bub3+bub4);
1382 return 16.*p1*p2*p3*p4/den2;
1397 libmesh_error_msg(
"Invalid i = " << i);
1404 libmesh_assert_less (i, 14);
1407 const Real zeta1 = p(0);
1408 const Real zeta2 = p(1);
1409 const Real zeta3 = p(2);
1410 const Real zeta0 = 1. - zeta1 - zeta2 - zeta3;
1413 const Real bubble_012 = zeta0*zeta1*zeta2;
1414 const Real bubble_013 = zeta0*zeta1*zeta3;
1415 const Real bubble_123 = zeta1*zeta2*zeta3;
1416 const Real bubble_023 = zeta0*zeta2*zeta3;
1421 return zeta0*(2.*zeta0 - 1.) + 3.*(bubble_012+bubble_013+bubble_023);
1424 return zeta1*(2.*zeta1 - 1.) + 3.*(bubble_012+bubble_013+bubble_123);
1427 return zeta2*(2.*zeta2 - 1.) + 3.*(bubble_012+bubble_023+bubble_123);
1430 return zeta3*(2.*zeta3 - 1.) + 3.*(bubble_013+bubble_023+bubble_123);
1433 return 4.*zeta0*zeta1 - 12.*(bubble_012+bubble_013);
1436 return 4.*zeta1*zeta2 - 12.*(bubble_012+bubble_123);
1439 return 4.*zeta2*zeta0 - 12.*(bubble_012+bubble_023);
1442 return 4.*zeta0*zeta3 - 12.*(bubble_013+bubble_023);
1445 return 4.*zeta1*zeta3 - 12.*(bubble_013+bubble_123);
1448 return 4.*zeta2*zeta3 - 12.*(bubble_023+bubble_123);
1451 return 27.*bubble_012;
1454 return 27.*bubble_013;
1457 return 27.*bubble_123;
1460 return 27.*bubble_023;
1463 libmesh_error_msg(
"Invalid i = " << i);
1474 libmesh_error_msg(
"ERROR: Unsupported 3D FE order!: " << order);
1477 #else // LIBMESH_DIM != 3 1479 libmesh_not_implemented();
1485 template <FEFamily T>
1489 const unsigned int i,
1490 const unsigned int j,
1493 #if LIBMESH_DIM == 3 1495 libmesh_assert_less (j, 3);
1509 libmesh_assert_less (i, 8);
1512 const Real xi = p(0);
1513 const Real eta = p(1);
1514 const Real zeta = p(2);
1516 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0};
1517 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1};
1518 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1};
1538 libmesh_error_msg(
"Invalid j = " << j);
1547 libmesh_assert_less (i, 4);
1550 const Real dzeta0dxi = -1.;
1551 const Real dzeta1dxi = 1.;
1552 const Real dzeta2dxi = 0.;
1553 const Real dzeta3dxi = 0.;
1555 const Real dzeta0deta = -1.;
1556 const Real dzeta1deta = 0.;
1557 const Real dzeta2deta = 1.;
1558 const Real dzeta3deta = 0.;
1560 const Real dzeta0dzeta = -1.;
1561 const Real dzeta1dzeta = 0.;
1562 const Real dzeta2dzeta = 0.;
1563 const Real dzeta3dzeta = 1.;
1585 libmesh_error_msg(
"Invalid i = " << i);
1607 libmesh_error_msg(
"Invalid i = " << i);
1629 libmesh_error_msg(
"Invalid i = " << i);
1634 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
1645 libmesh_assert_less (i, 6);
1650 Point p2d(p(0),p(1));
1654 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1};
1655 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2};
1675 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
1685 libmesh_assert_less (i, 5);
1687 const Real xi = p(0);
1688 const Real eta = p(1);
1689 const Real zeta = p(2);
1690 const Real eps = 1.e-35;
1699 return .25*(zeta + eta - 1.)/((1. - zeta) + eps);
1702 return -.25*(zeta + eta - 1.)/((1. - zeta) + eps);
1705 return -.25*(zeta - eta - 1.)/((1. - zeta) + eps);
1708 return .25*(zeta - eta - 1.)/((1. - zeta) + eps);
1714 libmesh_error_msg(
"Invalid i = " << i);
1723 return .25*(zeta + xi - 1.)/((1. - zeta) + eps);
1726 return .25*(zeta - xi - 1.)/((1. - zeta) + eps);
1729 return -.25*(zeta - xi - 1.)/((1. - zeta) + eps);
1732 return -.25*(zeta + xi - 1.)/((1. - zeta) + eps);
1738 libmesh_error_msg(
"Invalid i = " << i);
1748 num = zeta*(2. - zeta) - 1.,
1749 den = (1. - zeta + eps)*(1. - zeta + eps);
1755 return .25*(num + xi*eta)/den;
1759 return .25*(num - xi*eta)/den;
1765 libmesh_error_msg(
"Invalid i = " << i);
1770 libmesh_error_msg(
"Invalid j = " << j);
1778 libmesh_error_msg(
"Code (see stack trace) used an outdated FE function overload.\n" 1779 "Shape functions on a polyhedron are not defined by ElemType alone.");
1783 const C0Polyhedron & poly = *cast_ptr<const C0Polyhedron *>(elem);
1791 libmesh_assert_less(s, poly.n_subelements());
1793 const auto subtet = poly.subelement(s);
1797 Real du_da = 0, du_db = 0, du_dc = 0;
1800 const int nodei = i;
1801 if (nodei == subtet[0])
1802 du_da = du_db = du_dc = -1;
1803 else if (nodei == subtet[1])
1805 else if (nodei == subtet[2])
1807 else if (nodei == subtet[3])
1819 const auto master_points = poly.master_subelement(s);
1822 master_points[1](0) - master_points[0](0), master_points[2](0) - master_points[0](0), master_points[3](0) - master_points[0](0),
1823 master_points[1](1) - master_points[0](1), master_points[2](1) - master_points[0](1), master_points[3](1) - master_points[0](1),
1824 master_points[1](2) - master_points[0](2), master_points[2](2) - master_points[0](2), master_points[3](2) - master_points[0](2));
1830 const Real jac = dXi_dA.det();
1837 const Real da_dxi = (dXi_dA(2,2)*dXi_dA(1,1) - dXi_dA(2,1)*dXi_dA(1,2)) / jac;
1838 const Real db_dxi = -(dXi_dA(2,2)*dXi_dA(1,0) - dXi_dA(2,0)*dXi_dA(1,2)) / jac;
1839 const Real dc_dxi = (dXi_dA(2,1)*dXi_dA(1,0) - dXi_dA(2,0)*dXi_dA(1,1)) / jac;
1840 return du_da*da_dxi + du_db*db_dxi + du_dc*dc_dxi;
1845 const Real da_deta = -(dXi_dA(2,2)*dXi_dA(0,1) - dXi_dA(2,1)*dXi_dA(0,2) ) / jac;
1846 const Real db_deta = (dXi_dA(2,2)*dXi_dA(0,0) - dXi_dA(2,0)*dXi_dA(0,2)) / jac;
1847 const Real dc_deta = -(dXi_dA(2,1)*dXi_dA(0,0) - dXi_dA(2,0)*dXi_dA(0,1)) / jac;
1848 return du_da*da_deta + du_db*db_deta + du_dc*dc_deta;
1853 const Real da_dzeta = (dXi_dA(1,2)*dXi_dA(0,1) - dXi_dA(1,1)*dXi_dA(0,2) ) / jac;
1854 const Real db_dzeta = -(dXi_dA(1,2)*dXi_dA(0,0) - dXi_dA(1,0)*dXi_dA(0,2)) / jac;
1855 const Real dc_dzeta = (dXi_dA(1,1)*dXi_dA(0,0) - dXi_dA(1,0)*dXi_dA(0,1)) / jac;
1856 return du_da*da_dzeta + du_db*db_dzeta + du_dc*dc_dzeta;
1859 libmesh_error_msg(
"ERROR: Invalid derivative index j = " << j);
1878 libmesh_assert_less (i, 20);
1880 const Real xi = p(0);
1881 const Real eta = p(1);
1882 const Real zeta = p(2);
1886 const Real x = .5*(xi + 1.);
1887 const Real y = .5*(eta + 1.);
1888 const Real z = .5*(zeta + 1.);
1900 return .5*(1. - y)*(1. - z)*((1. - x)*(-2.) +
1901 (-1.)*(1. - 2.*x - 2.*y - 2.*z));
1904 return .5*(1. - y)*(1. - z)*(x*(2.) +
1905 (1.)*(2.*x - 2.*y - 2.*z - 1.));
1908 return .5*y*(1. - z)*(x*(2.) +
1909 (1.)*(2.*x + 2.*y - 2.*z - 3.));
1912 return .5*y*(1. - z)*((1. - x)*(-2.) +
1913 (-1.)*(2.*y - 2.*x - 2.*z - 1.));
1916 return .5*(1. - y)*z*((1. - x)*(-2.) +
1917 (-1.)*(2.*z - 2.*x - 2.*y - 1.));
1920 return .5*(1. - y)*z*(x*(2.) +
1921 (1.)*(2.*x - 2.*y + 2.*z - 3.));
1924 return .5*y*z*(x*(2.) +
1925 (1.)*(2.*x + 2.*y + 2.*z - 5.));
1928 return .5*y*z*((1. - x)*(-2.) +
1929 (-1.)*(2.*y - 2.*x + 2.*z - 3.));
1932 return 2.*(1. - y)*(1. - z)*(1. - 2.*x);
1935 return 2.*y*(1. - y)*(1. - z);
1938 return 2.*y*(1. - z)*(1. - 2.*x);
1941 return 2.*y*(1. - y)*(1. - z)*(-1.);
1944 return 2.*(1. - y)*z*(1. - z)*(-1.);
1947 return 2.*(1. - y)*z*(1. - z);
1950 return 2.*y*z*(1. - z);
1953 return 2.*y*z*(1. - z)*(-1.);
1956 return 2.*(1. - y)*z*(1. - 2.*x);
1959 return 2.*y*(1. - y)*z;
1962 return 2.*y*z*(1. - 2.*x);
1965 return 2.*y*(1. - y)*z*(-1.);
1968 libmesh_error_msg(
"Invalid i = " << i);
1977 return .5*(1. - x)*(1. - z)*((1. - y)*(-2.) +
1978 (-1.)*(1. - 2.*x - 2.*y - 2.*z));
1981 return .5*x*(1. - z)*((1. - y)*(-2.) +
1982 (-1.)*(2.*x - 2.*y - 2.*z - 1.));
1985 return .5*x*(1. - z)*(y*(2.) +
1986 (1.)*(2.*x + 2.*y - 2.*z - 3.));
1989 return .5*(1. - x)*(1. - z)*(y*(2.) +
1990 (1.)*(2.*y - 2.*x - 2.*z - 1.));
1993 return .5*(1. - x)*z*((1. - y)*(-2.) +
1994 (-1.)*(2.*z - 2.*x - 2.*y - 1.));
1997 return .5*x*z*((1. - y)*(-2.) +
1998 (-1.)*(2.*x - 2.*y + 2.*z - 3.));
2001 return .5*x*z*(y*(2.) +
2002 (1.)*(2.*x + 2.*y + 2.*z - 5.));
2005 return .5*(1. - x)*z*(y*(2.) +
2006 (1.)*(2.*y - 2.*x + 2.*z - 3.));
2009 return 2.*x*(1. - x)*(1. - z)*(-1.);
2012 return 2.*x*(1. - z)*(1. - 2.*y);
2015 return 2.*x*(1. - x)*(1. - z);
2018 return 2.*(1. - x)*(1. - z)*(1. - 2.*y);
2021 return 2.*(1. - x)*z*(1. - z)*(-1.);
2024 return 2.*x*z*(1. - z)*(-1.);
2027 return 2.*x*z*(1. - z);
2030 return 2.*(1. - x)*z*(1. - z);
2033 return 2.*x*(1. - x)*z*(-1.);
2036 return 2.*x*z*(1. - 2.*y);
2039 return 2.*x*(1. - x)*z;
2042 return 2.*(1. - x)*z*(1. - 2.*y);
2045 libmesh_error_msg(
"Invalid i = " << i);
2054 return .5*(1. - x)*(1. - y)*((1. - z)*(-2.) +
2055 (-1.)*(1. - 2.*x - 2.*y - 2.*z));
2058 return .5*x*(1. - y)*((1. - z)*(-2.) +
2059 (-1.)*(2.*x - 2.*y - 2.*z - 1.));
2062 return .5*x*y*((1. - z)*(-2.) +
2063 (-1.)*(2.*x + 2.*y - 2.*z - 3.));
2066 return .5*(1. - x)*y*((1. - z)*(-2.) +
2067 (-1.)*(2.*y - 2.*x - 2.*z - 1.));
2070 return .5*(1. - x)*(1. - y)*(z*(2.) +
2071 (1.)*(2.*z - 2.*x - 2.*y - 1.));
2074 return .5*x*(1. - y)*(z*(2.) +
2075 (1.)*(2.*x - 2.*y + 2.*z - 3.));
2078 return .5*x*y*(z*(2.) +
2079 (1.)*(2.*x + 2.*y + 2.*z - 5.));
2082 return .5*(1. - x)*y*(z*(2.) +
2083 (1.)*(2.*y - 2.*x + 2.*z - 3.));
2086 return 2.*x*(1. - x)*(1. - y)*(-1.);
2089 return 2.*x*y*(1. - y)*(-1.);
2092 return 2.*x*(1. - x)*y*(-1.);
2095 return 2.*(1. - x)*y*(1. - y)*(-1.);
2098 return 2.*(1. - x)*(1. - y)*(1. - 2.*z);
2101 return 2.*x*(1. - y)*(1. - 2.*z);
2104 return 2.*x*y*(1. - 2.*z);
2107 return 2.*(1. - x)*y*(1. - 2.*z);
2110 return 2.*x*(1. - x)*(1. - y);
2113 return 2.*x*y*(1. - y);
2116 return 2.*x*(1. - x)*y;
2119 return 2.*(1. - x)*y*(1. - y);
2122 libmesh_error_msg(
"Invalid i = " << i);
2126 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
2133 "High order on first order elements only supported for L2 families");
2134 libmesh_fallthrough();
2137 libmesh_assert_less (i, 27);
2140 const Real xi = p(0);
2141 const Real eta = p(1);
2142 const Real zeta = p(2);
2148 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 0, 2, 2, 1, 2, 0, 2, 2};
2149 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 2, 0, 2, 1, 2, 2, 2};
2150 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1, 0, 0, 0, 0, 2, 2, 2, 2, 1, 1, 1, 1, 0, 2, 2, 2, 2, 1, 2};
2170 libmesh_error_msg(
"Invalid j = " << j);
2177 "High order on first order elements only supported for L2 families");
2178 libmesh_fallthrough();
2182 libmesh_assert_less (i, 10);
2185 const Real zeta1 = p(0);
2186 const Real zeta2 = p(1);
2187 const Real zeta3 = p(2);
2188 const Real zeta0 = 1. - zeta1 - zeta2 - zeta3;
2190 const Real dzeta0dxi = -1.;
2191 const Real dzeta1dxi = 1.;
2192 const Real dzeta2dxi = 0.;
2193 const Real dzeta3dxi = 0.;
2195 const Real dzeta0deta = -1.;
2196 const Real dzeta1deta = 0.;
2197 const Real dzeta2deta = 1.;
2198 const Real dzeta3deta = 0.;
2200 const Real dzeta0dzeta = -1.;
2201 const Real dzeta1dzeta = 0.;
2202 const Real dzeta2dzeta = 0.;
2203 const Real dzeta3dzeta = 1.;
2213 return (4.*zeta0 - 1.)*dzeta0dxi;
2216 return (4.*zeta1 - 1.)*dzeta1dxi;
2219 return (4.*zeta2 - 1.)*dzeta2dxi;
2222 return (4.*zeta3 - 1.)*dzeta3dxi;
2225 return 4.*(zeta0*dzeta1dxi + dzeta0dxi*zeta1);
2228 return 4.*(zeta1*dzeta2dxi + dzeta1dxi*zeta2);
2231 return 4.*(zeta0*dzeta2dxi + dzeta0dxi*zeta2);
2234 return 4.*(zeta0*dzeta3dxi + dzeta0dxi*zeta3);
2237 return 4.*(zeta1*dzeta3dxi + dzeta1dxi*zeta3);
2240 return 4.*(zeta2*dzeta3dxi + dzeta2dxi*zeta3);
2243 libmesh_error_msg(
"Invalid i = " << i);
2253 return (4.*zeta0 - 1.)*dzeta0deta;
2256 return (4.*zeta1 - 1.)*dzeta1deta;
2259 return (4.*zeta2 - 1.)*dzeta2deta;
2262 return (4.*zeta3 - 1.)*dzeta3deta;
2265 return 4.*(zeta0*dzeta1deta + dzeta0deta*zeta1);
2268 return 4.*(zeta1*dzeta2deta + dzeta1deta*zeta2);
2271 return 4.*(zeta0*dzeta2deta + dzeta0deta*zeta2);
2274 return 4.*(zeta0*dzeta3deta + dzeta0deta*zeta3);
2277 return 4.*(zeta1*dzeta3deta + dzeta1deta*zeta3);
2280 return 4.*(zeta2*dzeta3deta + dzeta2deta*zeta3);
2283 libmesh_error_msg(
"Invalid i = " << i);
2293 return (4.*zeta0 - 1.)*dzeta0dzeta;
2296 return (4.*zeta1 - 1.)*dzeta1dzeta;
2299 return (4.*zeta2 - 1.)*dzeta2dzeta;
2302 return (4.*zeta3 - 1.)*dzeta3dzeta;
2305 return 4.*(zeta0*dzeta1dzeta + dzeta0dzeta*zeta1);
2308 return 4.*(zeta1*dzeta2dzeta + dzeta1dzeta*zeta2);
2311 return 4.*(zeta0*dzeta2dzeta + dzeta0dzeta*zeta2);
2314 return 4.*(zeta0*dzeta3dzeta + dzeta0dzeta*zeta3);
2317 return 4.*(zeta1*dzeta3dzeta + dzeta1dzeta*zeta3);
2320 return 4.*(zeta2*dzeta3dzeta + dzeta2dzeta*zeta3);
2323 libmesh_error_msg(
"Invalid i = " << i);
2328 libmesh_error_msg(
"Invalid j = " << j);
2336 libmesh_assert_less (i, 15);
2338 const Real xi = p(0);
2339 const Real eta = p(1);
2340 const Real zeta = p(2);
2350 return (2.*xi + 2.*eta + 0.5*zeta - 1.)*(1. - zeta);
2352 return (2.*xi - 1. - 0.5*zeta)*(1. - zeta);
2356 return (2.*xi + 2.*eta - 0.5*zeta - 1.)*(1. + zeta);
2358 return (2.*xi - 1. + 0.5*zeta)*(1. + zeta);
2362 return (4.*xi + 2.*eta - 2.)*(zeta - 1.);
2364 return -2.*(zeta - 1.)*eta;
2366 return 2.*(zeta - 1.)*eta;
2368 return (zeta - 1.)*(1. + zeta);
2370 return (1. - zeta)*(1. + zeta);
2374 return (-4.*xi - 2.*eta + 2.)*(1. + zeta);
2376 return 2.*(1. + zeta)*eta;
2378 return -2.*(1. + zeta)*eta;
2380 libmesh_error_msg(
"Invalid i = " << i);
2390 return (2.*xi + 2.*eta + 0.5*zeta - 1.)*(1. - zeta);
2394 return (2.*eta - 1. - 0.5*zeta)*(1. - zeta);
2396 return (2.*xi + 2.*eta - 0.5*zeta - 1.)*(1. + zeta);
2400 return (2.*eta - 1. + 0.5*zeta)*(1. + zeta);
2402 return 2.*(zeta - 1.)*xi;
2404 return 2.*(1. - zeta)*xi;
2406 return (2.*xi + 4.*eta - 2.)*(zeta - 1.);
2408 return (zeta - 1.)*(1. + zeta);
2412 return (1. - zeta)*(1. + zeta);
2414 return -2.*(1. + zeta)*xi;
2416 return 2.*(1. + zeta)*xi;
2418 return (-2.*xi - 4.*eta + 2.)*(1. + zeta);
2421 libmesh_error_msg(
"Invalid i = " << i);
2431 return (-xi - eta - zeta + 0.5)*(xi + eta - 1.);
2433 return -0.5*xi*(2.*xi - 1. - 2.*zeta);
2435 return -0.5*eta*(2.*eta - 1. - 2.*zeta);
2437 return (xi + eta - zeta - 0.5)*(xi + eta - 1.);
2439 return 0.5*xi*(2.*xi - 1. + 2.*zeta);
2441 return 0.5*eta*(2.*eta - 1. + 2.*zeta);
2443 return 2.*xi*(xi + eta - 1.);
2447 return 2.*eta*(xi + eta - 1.);
2449 return 2.*zeta*(xi + eta - 1.);
2453 return -2.*eta*zeta;
2455 return 2.*xi*(1. - xi - eta);
2459 return 2.*eta*(1. - xi - eta);
2462 libmesh_error_msg(
"Invalid i = " << i);
2467 libmesh_error_msg(
"Invalid j = " << j);
2476 "High order on first order elements only supported for L2 families");
2477 libmesh_fallthrough();
2482 libmesh_assert_less (i, 18);
2487 Point p2d(p(0),p(1));
2491 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2};
2492 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2, 3, 4, 5, 0, 1, 2, 3, 4, 5, 3, 4, 5};
2512 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
2521 libmesh_assert_less (i, 13);
2523 const Real xi = p(0);
2524 const Real eta = p(1);
2525 const Real zeta = p(2);
2526 const Real eps = 1.e-35;
2531 den = (-1. + zeta + eps),
2545 return 0.25*(-zeta - eta + 2.*eta*zeta - 2.*xi + 2.*zeta*xi + 2.*eta*xi + zeta2 + eta2)/den;
2548 return -0.25*(-zeta - eta + 2.*eta*zeta + 2.*xi - 2.*zeta*xi - 2.*eta*xi + zeta2 + eta2)/den;
2551 return -0.25*(-zeta + eta - 2.*eta*zeta + 2.*xi - 2.*zeta*xi + 2.*eta*xi + zeta2 + eta2)/den;
2554 return 0.25*(-zeta + eta - 2.*eta*zeta - 2.*xi + 2.*zeta*xi - 2.*eta*xi + zeta2 + eta2)/den;
2560 return -(-1. + eta + zeta)*xi/den;
2563 return 0.5*(-1. + eta + zeta)*(1. + eta - zeta)/den;
2566 return (1. + eta - zeta)*xi/den;
2569 return -0.5*(-1. + eta + zeta)*(1. + eta - zeta)/den;
2572 return -(-1. + eta + zeta)*zeta/den;
2575 return (-1. + eta + zeta)*zeta/den;
2578 return -(1. + eta - zeta)*zeta/den;
2581 return (1. + eta - zeta)*zeta/den;
2584 libmesh_error_msg(
"Invalid i = " << i);
2592 return 0.25*(-zeta - 2.*eta + 2.*eta*zeta - xi + 2.*zeta*xi + 2.*eta*xi + zeta2 + xi2)/den;
2595 return -0.25*(zeta + 2.*eta - 2.*eta*zeta - xi + 2.*zeta*xi + 2.*eta*xi - zeta2 - xi2)/den;
2598 return -0.25*(-zeta + 2.*eta - 2.*eta*zeta + xi - 2.*zeta*xi + 2.*eta*xi + zeta2 + xi2)/den;
2601 return 0.25*(zeta - 2.*eta + 2.*eta*zeta + xi - 2.*zeta*xi + 2.*eta*xi - zeta2 - xi2)/den;
2607 return -0.5*(-1. + xi + zeta)*(1. + xi - zeta)/den;
2610 return (1. + xi - zeta)*eta/den;
2613 return 0.5*(-1. + xi + zeta)*(1. + xi - zeta)/den;
2616 return -(-1. + xi + zeta)*eta/den;
2619 return -(-1. + xi + zeta)*zeta/den;
2622 return (1. + xi - zeta)*zeta/den;
2625 return -(1. + xi - zeta)*zeta/den;
2628 return (-1. + xi + zeta)*zeta/den;
2631 libmesh_error_msg(
"Invalid i = " << i);
2640 return -0.25*(xi + eta + 1.)*(-1. + 2.*zeta - zeta2 + eta*xi)/den2;
2643 return 0.25*(eta - xi + 1.)*(1. - 2.*zeta + zeta2 + eta*xi)/den2;
2646 return 0.25*(xi + eta - 1.)*(-1. + 2.*zeta - zeta2 + eta*xi)/den2;
2649 return -0.25*(eta - xi - 1.)*(1. - 2.*zeta + zeta2 + eta*xi)/den2;
2652 return 4.*zeta - 1.;
2655 return 0.5*(-2 + eta + 6.*zeta + eta*xi2 + eta*zeta2 - 6.*zeta2 + 2.*zeta3 - 2.*eta*zeta)/den2;
2658 return -0.5*(2 - 6.*zeta + xi + xi*zeta2 + eta2*xi + 6.*zeta2 - 2.*zeta3 - 2.*zeta*xi)/den2;
2661 return -0.5*(2 + eta - 6.*zeta + eta*xi2 + eta*zeta2 + 6.*zeta2 - 2.*zeta3 - 2.*eta*zeta)/den2;
2664 return 0.5*(-2 + 6.*zeta + xi + xi*zeta2 + eta2*xi - 6.*zeta2 + 2.*zeta3 - 2.*zeta*xi)/den2;
2667 return (1. - eta - 4.*zeta - xi - xi*zeta2 - eta*zeta2 + eta*xi + 5.*zeta2 - 2.*zeta3 + 2.*eta*zeta + 2.*zeta*xi)/den2;
2670 return -(-1. + eta + 4.*zeta - xi - xi*zeta2 + eta*zeta2 + eta*xi - 5.*zeta2 + 2.*zeta3 - 2.*eta*zeta + 2.*zeta*xi)/den2;
2673 return (1. + eta - 4.*zeta + xi + xi*zeta2 + eta*zeta2 + eta*xi + 5.*zeta2 - 2.*zeta3 - 2.*eta*zeta - 2.*zeta*xi)/den2;
2676 return -(-1. - eta + 4.*zeta + xi + xi*zeta2 - eta*zeta2 + eta*xi - 5.*zeta2 + 2.*zeta3 + 2.*eta*zeta - 2.*zeta*xi)/den2;
2679 libmesh_error_msg(
"Invalid i = " << i);
2684 libmesh_error_msg(
"Invalid j = " << j);
2693 "High order on first order elements only supported for L2 families");
2694 libmesh_fallthrough();
2698 libmesh_assert_less (i, 14);
2700 const Real xi = p(0);
2701 const Real eta = p(1);
2702 const Real zeta = p(2);
2703 const Real eps = 1.e-35;
2708 p1 = 0.5*(1. - eta - zeta),
2709 p2 = 0.5*(1. + xi - zeta),
2710 p3 = 0.5*(1. + eta - zeta),
2711 p4 = 0.5*(1. - xi - zeta);
2716 den = (-1. + zeta + eps),
2727 return 0.5*p1*(-xi*eta + zeta - zeta*zeta + 2.*p4*eta)/den2;
2730 return -0.5*p1*(xi*eta + zeta - zeta*zeta + 2.*p2*eta)/den2;
2733 return 0.5*p3*(xi*eta - zeta + zeta*zeta + 2.*p2*eta)/den2;
2736 return -0.5*p3*(-xi*eta - zeta + zeta*zeta + 2.*p4*eta)/den2;
2742 return 2.*p1*eta*xi/den2;
2745 return 2.*p1*p3*(xi + 2.*p2)/den2;
2748 return -2.*p3*eta*xi/den2;
2751 return -2.*p1*p3*(-xi + 2.*p4)/den2;
2754 return 2.*p1*zeta/den;
2757 return -2.*p1*zeta/den;
2760 return -2.*p3*zeta/den;
2763 return 2.*p3*zeta/den;
2766 return -8.*p1*p3*xi/den2;
2769 libmesh_error_msg(
"Invalid i = " << i);
2777 return -0.5*p4*(xi*eta - zeta + zeta*zeta - 2.*p1*xi)/den2;
2780 return 0.5*p2*(xi*eta + zeta - zeta*zeta - 2.*p1*xi)/den2;
2783 return 0.5*p2*(xi*eta - zeta + zeta*zeta + 2.*p3*xi)/den2;
2786 return -0.5*p4*(xi*eta + zeta - zeta*zeta + 2.*p3*xi)/den2;
2792 return 2.*p2*p4*(eta - 2.*p1)/den2;
2795 return -2.*p2*xi*eta/den2;
2798 return 2.*p2*p4*(eta + 2.*p3)/den2;
2801 return 2.*p4*xi*eta/den2;
2804 return 2.*p4*zeta/den;
2807 return 2.*p2*zeta/den;
2810 return -2.*p2*zeta/den;
2813 return -2.*p4*zeta/den;
2816 return -8.*p2*p4*eta/den2;
2819 libmesh_error_msg(
"Invalid i = " << i);
2829 return -0.5*p1*(xi*eta - zeta + zeta*zeta)/den2
2830 - 0.5*p4*(xi*eta - zeta + zeta*zeta)/den2
2831 + p4*p1*(2.*zeta - 1)/den2
2832 - 2.*p4*p1*(xi*eta - zeta + zeta*zeta)/den3;
2835 return 0.5*p2*(xi*eta + zeta - zeta*zeta)/den2
2836 + 0.5*p1*(xi*eta + zeta - zeta*zeta)/den2
2837 - p1*p2*(1 - 2.*zeta)/den2
2838 + 2.*p1*p2*(xi*eta + zeta - zeta*zeta)/den3;
2841 return -0.5*p3*(xi*eta - zeta + zeta*zeta)/den2
2842 - 0.5*p2*(xi*eta - zeta + zeta*zeta)/den2
2843 + p2*p3*(2.*zeta - 1)/den2
2844 - 2.*p2*p3*(xi*eta - zeta + zeta*zeta)/den3;
2847 return 0.5*p4*(xi*eta + zeta - zeta*zeta)/den2
2848 + 0.5*p3*(xi*eta + zeta - zeta*zeta)/den2
2849 - p3*p4*(1 - 2.*zeta)/den2
2850 + 2.*p3*p4*(xi*eta + zeta - zeta*zeta)/den3;
2853 return 4.*zeta - 1.;
2856 return 2.*p4*p1*eta/den2
2859 + 8.*p2*p1*p4*eta/den3;
2862 return -2.*p2*p3*xi/den2
2865 - 8.*p1*p2*p3*xi/den3;
2868 return -2.*p3*p4*eta/den2
2871 - 8.*p2*p3*p4*eta/den3;
2874 return 2.*p4*p1*xi/den2
2877 + 8.*p3*p4*p1*xi/den3;
2880 return 2.*p4*zeta/den
2883 + 4.*p1*p4*zeta/den2;
2886 return 2.*p1*zeta/den
2889 + 4.*p2*p1*zeta/den2;
2892 return 2.*p2*zeta/den
2895 + 4.*p3*p2*zeta/den2;
2898 return 2.*p3*zeta/den
2901 + 4.*p4*p3*zeta/den2;
2904 return -8.*p2*p3*p4/den2
2908 - 32.*p1*p2*p3*p4/den3;
2911 libmesh_error_msg(
"Invalid i = " << i);
2916 libmesh_error_msg(
"Invalid j = " << j);
2933 libmesh_assert_less (i, 14);
2936 const Real zeta1 = p(0);
2937 const Real zeta2 = p(1);
2938 const Real zeta3 = p(2);
2939 const Real zeta0 = 1. - zeta1 - zeta2 - zeta3;
2941 const Real dzeta0dxi = -1.;
2942 const Real dzeta1dxi = 1.;
2943 const Real dzeta2dxi = 0.;
2944 const Real dzeta3dxi = 0.;
2945 const Real dbubble012dxi = (zeta0-zeta1)*zeta2;
2946 const Real dbubble013dxi = (zeta0-zeta1)*zeta3;
2947 const Real dbubble123dxi = zeta2*zeta3;
2948 const Real dbubble023dxi = -zeta2*zeta3;
2950 const Real dzeta0deta = -1.;
2951 const Real dzeta1deta = 0.;
2952 const Real dzeta2deta = 1.;
2953 const Real dzeta3deta = 0.;
2954 const Real dbubble012deta = (zeta0-zeta2)*zeta1;
2955 const Real dbubble013deta = -zeta1*zeta3;
2956 const Real dbubble123deta = zeta1*zeta3;
2957 const Real dbubble023deta = (zeta0-zeta2)*zeta3;
2959 const Real dzeta0dzeta = -1.;
2960 const Real dzeta1dzeta = 0.;
2961 const Real dzeta2dzeta = 0.;
2962 const Real dzeta3dzeta = 1.;
2963 const Real dbubble012dzeta = -zeta1*zeta2;
2964 const Real dbubble013dzeta = (zeta0-zeta3)*zeta1;
2965 const Real dbubble123dzeta = zeta1*zeta2;
2966 const Real dbubble023dzeta = (zeta0-zeta3)*zeta2;
2976 return (4.*zeta0 - 1.)*dzeta0dxi + 3.*(dbubble012dxi+dbubble013dxi+dbubble023dxi);
2979 return (4.*zeta1 - 1.)*dzeta1dxi + 3.*(dbubble012dxi+dbubble013dxi+dbubble123dxi);
2982 return (4.*zeta2 - 1.)*dzeta2dxi + 3.*(dbubble012dxi+dbubble023dxi+dbubble123dxi);
2985 return (4.*zeta3 - 1.)*dzeta3dxi + 3.*(dbubble013dxi+dbubble023dxi+dbubble123dxi);
2988 return 4.*(zeta0*dzeta1dxi + dzeta0dxi*zeta1) - 12.*(dbubble012dxi+dbubble013dxi);
2991 return 4.*(zeta1*dzeta2dxi + dzeta1dxi*zeta2) - 12.*(dbubble012dxi+dbubble123dxi);
2994 return 4.*(zeta0*dzeta2dxi + dzeta0dxi*zeta2) - 12.*(dbubble012dxi+dbubble023dxi);
2997 return 4.*(zeta0*dzeta3dxi + dzeta0dxi*zeta3) - 12.*(dbubble013dxi+dbubble023dxi);
3000 return 4.*(zeta1*dzeta3dxi + dzeta1dxi*zeta3) - 12.*(dbubble013dxi+dbubble123dxi);
3003 return 4.*(zeta2*dzeta3dxi + dzeta2dxi*zeta3) - 12.*(dbubble023dxi+dbubble123dxi);
3006 return 27.*dbubble012dxi;
3009 return 27.*dbubble013dxi;
3012 return 27.*dbubble123dxi;
3015 return 27.*dbubble023dxi;
3018 libmesh_error_msg(
"Invalid i = " << i);
3028 return (4.*zeta0 - 1.)*dzeta0deta + 3.*(dbubble012deta+dbubble013deta+dbubble023deta);;
3031 return (4.*zeta1 - 1.)*dzeta1deta + 3.*(dbubble012deta+dbubble013deta+dbubble123deta);
3034 return (4.*zeta2 - 1.)*dzeta2deta + 3.*(dbubble012deta+dbubble023deta+dbubble123deta);
3037 return (4.*zeta3 - 1.)*dzeta3deta + 3.*(dbubble013deta+dbubble023deta+dbubble123deta);
3040 return 4.*(zeta0*dzeta1deta + dzeta0deta*zeta1) - 12.*(dbubble012deta+dbubble013deta);
3043 return 4.*(zeta1*dzeta2deta + dzeta1deta*zeta2) - 12.*(dbubble012deta+dbubble123deta);
3046 return 4.*(zeta0*dzeta2deta + dzeta0deta*zeta2) - 12.*(dbubble012deta+dbubble023deta);
3049 return 4.*(zeta0*dzeta3deta + dzeta0deta*zeta3) - 12.*(dbubble013deta+dbubble023deta);
3052 return 4.*(zeta1*dzeta3deta + dzeta1deta*zeta3) - 12.*(dbubble013deta+dbubble123deta);
3055 return 4.*(zeta2*dzeta3deta + dzeta2deta*zeta3) - 12.*(dbubble023deta+dbubble123deta);
3058 return 27.*dbubble012deta;
3061 return 27.*dbubble013deta;
3064 return 27.*dbubble123deta;
3067 return 27.*dbubble023deta;
3070 libmesh_error_msg(
"Invalid i = " << i);
3080 return (4.*zeta0 - 1.)*dzeta0dzeta + 3.*(dbubble012dzeta+dbubble013dzeta+dbubble023dzeta);
3083 return (4.*zeta1 - 1.)*dzeta1dzeta + 3.*(dbubble012dzeta+dbubble013dzeta+dbubble123dzeta);
3086 return (4.*zeta2 - 1.)*dzeta2dzeta + 3.*(dbubble012dzeta+dbubble023dzeta+dbubble123dzeta);
3089 return (4.*zeta3 - 1.)*dzeta3dzeta + 3.*(dbubble013dzeta+dbubble023dzeta+dbubble123dzeta);
3092 return 4.*(zeta0*dzeta1dzeta + dzeta0dzeta*zeta1) - 12.*(dbubble012dzeta+dbubble013dzeta);
3095 return 4.*(zeta1*dzeta2dzeta + dzeta1dzeta*zeta2) - 12.*(dbubble012dzeta+dbubble123dzeta);
3098 return 4.*(zeta0*dzeta2dzeta + dzeta0dzeta*zeta2) - 12.*(dbubble012dzeta+dbubble023dzeta);
3101 return 4.*(zeta0*dzeta3dzeta + dzeta0dzeta*zeta3) - 12.*(dbubble013dzeta+dbubble023dzeta);
3104 return 4.*(zeta1*dzeta3dzeta + dzeta1dzeta*zeta3) - 12.*(dbubble013dzeta+dbubble123dzeta);
3107 return 4.*(zeta2*dzeta3dzeta + dzeta2dzeta*zeta3) - 12.*(dbubble023dzeta+dbubble123dzeta);
3110 return 27.*dbubble012dzeta;
3113 return 27.*dbubble013dzeta;
3116 return 27.*dbubble123dzeta;
3119 return 27.*dbubble023dzeta;
3122 libmesh_error_msg(
"Invalid i = " << i);
3127 libmesh_error_msg(
"Invalid j = " << j);
3133 libmesh_assert_less (i, 20);
3136 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2, 0, 1, 2};
3138 Point p2d(p(0),p(1));
3172 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
3176 return mainval - bubbleval / 9;
3178 return mainval + bubbleval * (
Real(4) / 9);
3183 libmesh_assert_less (i, 21);
3188 Point p2d(p(0),p(1));
3192 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2, 0, 1, 2};
3193 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2, 3, 4, 5, 0, 1, 2, 3, 4, 5, 3, 4, 5, 6, 6, 6};
3213 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
3219 return fe_fdm_deriv(type, order, elem, i, j, p, fe_lagrange_3D_shape<T>);
3229 libmesh_error_msg(
"ERROR: Unsupported 3D FE order!: " << order);
3232 #else // LIBMESH_DIM != 3 3234 libmesh_not_implemented();
3240 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 3242 template <FEFamily T>
3243 Real fe_lagrange_3D_shape_second_deriv(
const ElemType type,
3246 const unsigned int i,
3247 const unsigned int j,
3250 #if LIBMESH_DIM == 3 3252 libmesh_assert_less (j, 6);
3279 libmesh_assert_less (i, 6);
3284 Point p2d(p(0),p(1));
3288 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1};
3289 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2};
3311 libmesh_error_msg(
"Invalid j = " << j);
3320 libmesh_assert_less (i, 5);
3322 const Real xi = p(0);
3323 const Real eta = p(1);
3324 const Real zeta = p(2);
3325 const Real eps = 1.e-35;
3340 return 0.25/(1. - zeta + eps);
3343 return -0.25/(1. - zeta + eps);
3347 libmesh_error_msg(
"Invalid i = " << i);
3353 Real den = (1. - zeta + eps)*(1. - zeta + eps);
3359 return 0.25*eta/den;
3362 return -0.25*eta/den;
3366 libmesh_error_msg(
"Invalid i = " << i);
3372 Real den = (1. - zeta + eps)*(1. - zeta + eps);
3381 return -0.25*xi/den;
3385 libmesh_error_msg(
"Invalid i = " << i);
3391 Real den = (1. - zeta + eps)*(1. - zeta + eps)*(1. - zeta + eps);
3397 return 0.5*xi*eta/den;
3400 return -0.5*xi*eta/den;
3404 libmesh_error_msg(
"Invalid i = " << i);
3409 libmesh_error_msg(
"Invalid j = " << j);
3418 libmesh_assert_less (i, 8);
3421 const Real xi = p(0);
3422 const Real eta = p(1);
3423 const Real zeta = p(2);
3425 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0};
3426 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1};
3427 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1};
3455 libmesh_error_msg(
"Invalid j = " << j);
3482 libmesh_assert_less (i, 20);
3484 const Real xi = p(0);
3485 const Real eta = p(1);
3486 const Real zeta = p(2);
3490 const Real x = .5*(xi + 1.);
3491 const Real y = .5*(eta + 1.);
3492 const Real z = .5*(zeta + 1.);
3502 return (1. - y) * (1. - z);
3505 return y * (1. - z);
3508 return (1. - y) * z;
3513 return -2. * (1. - y) * (1. - z);
3515 return -2. * y * (1. - z);
3517 return -2. * (1. - y) * z;
3530 libmesh_error_msg(
"Invalid i = " << i);
3538 return (1.25 - x - y - .5*z) * (1. - z);
3540 return (-x + y + .5*z - .25) * (1. - z);
3542 return (x + y - .5*z - .75) * (1. - z);
3544 return (-y + x + .5*z - .25) * (1. - z);
3546 return -.25*z * (4.*x + 4.*y - 2.*z - 3);
3548 return -.25*z * (-4.*y + 4.*x + 2.*z - 1.);
3550 return .25*z * (-5 + 4.*x + 4.*y + 2.*z);
3552 return .25*z * (4.*x - 4.*y - 2.*z + 1.);
3554 return (-1. + 2.*x) * (1. - z);
3556 return (1. - 2.*y) * (1. - z);
3558 return (1. - 2.*x) * (1. - z);
3560 return (-1. + 2.*y) * (1. - z);
3562 return z * (1. - z);
3564 return -z * (1. - z);
3566 return z * (1. - z);
3568 return -z * (1. - z);
3570 return (-1. + 2.*x) * z;
3572 return (1. - 2.*y) * z;
3574 return (1. - 2.*x) * z;
3576 return (-1. + 2.*y) * z;
3578 libmesh_error_msg(
"Invalid i = " << i);
3586 return (1. - x) * (1. - z);
3589 return x * (1. - z);
3592 return (1. - x) * z;
3597 return -2. * x * (1. - z);
3599 return -2. * (1. - x) * (1. - z);
3603 return -2. * (1. - x) * z;
3614 libmesh_error_msg(
"Invalid i = " << i);
3620 return (1.25 - x - .5*y - z) * (1. - y);
3622 return (-x + .5*y + z - .25) * (1. - y);
3624 return -.25*y * (2.*y + 4.*x - 4.*z - 1.);
3626 return -.25*y * (-2.*y + 4.*x + 4.*z - 3);
3628 return (-z + x + .5*y - .25) * (1. - y);
3630 return (x - .5*y + z - .75) * (1. - y);
3632 return .25*y * (2.*y + 4.*x + 4.*z - 5);
3634 return .25*y * (-2.*y + 4.*x - 4.*z + 1.);
3636 return (-1. + 2.*x) * (1. - y);
3638 return -y * (1. - y);
3640 return (-1. + 2.*x) * y;
3642 return y * (1. - y);
3644 return (-1. + 2.*z) * (1. - y);
3646 return (1. - 2.*z) * (1. - y);
3648 return (1. - 2.*z) * y;
3650 return (-1. + 2.*z) * y;
3652 return (1. - 2.*x) * (1. - y);
3654 return y * (1. - y);
3656 return (1. - 2.*x) * y;
3658 return -y * (1. - y);
3660 libmesh_error_msg(
"Invalid i = " << i);
3666 return (1.25 - .5*x - y - z) * (1. - x);
3668 return .25*x * (2.*x - 4.*y - 4.*z + 3.);
3670 return -.25*x * (2.*x + 4.*y - 4.*z - 1.);
3672 return (-y + .5*x + z - .25) * (1. - x);
3674 return (-z + .5*x + y - .25) * (1. - x);
3676 return -.25*x * (2.*x - 4.*y + 4.*z - 1.);
3678 return .25*x * (2.*x + 4.*y + 4.*z - 5);
3680 return (y - .5*x + z - .75) * (1. - x);
3682 return x * (1. - x);
3684 return (-1. + 2.*y) * x;
3686 return -x * (1. - x);
3688 return (-1. + 2.*y) * (1. - x);
3690 return (-1. + 2.*z) * (1. - x);
3692 return (-1. + 2.*z) * x;
3694 return (1. - 2.*z) * x;
3696 return (1. - 2.*z) * (1. - x);
3698 return -x * (1. - x);
3700 return (1. - 2.*y) * x;
3702 return x * (1. - x);
3704 return (1. - 2.*y) * (1. - x);
3706 libmesh_error_msg(
"Invalid i = " << i);
3713 return (1. - x) * (1. - y);
3716 return x * (1. - y);
3722 return (1. - x) * y;
3724 return -2. * (1. - x) * (1. - y);
3726 return -2. * x * (1. - y);
3730 return -2. * (1. - x) * y;
3741 libmesh_error_msg(
"Invalid i = " << i);
3744 libmesh_error_msg(
"Invalid j = " << j);
3751 "High order on first order elements only supported for L2 families");
3752 libmesh_fallthrough();
3755 libmesh_assert_less (i, 27);
3758 const Real xi = p(0);
3759 const Real eta = p(1);
3760 const Real zeta = p(2);
3766 static const unsigned int i0[] = {0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 0, 2, 2, 1, 2, 0, 2, 2};
3767 static const unsigned int i1[] = {0, 0, 1, 1, 0, 0, 1, 1, 0, 2, 1, 2, 0, 0, 1, 1, 0, 2, 1, 2, 2, 0, 2, 1, 2, 2, 2};
3768 static const unsigned int i2[] = {0, 0, 0, 0, 1, 1, 1, 1, 0, 0, 0, 0, 2, 2, 2, 2, 1, 1, 1, 1, 0, 2, 2, 2, 2, 1, 2};
3809 libmesh_error_msg(
"Invalid j = " << j);
3816 "High order on first order elements only supported for L2 families");
3817 libmesh_fallthrough();
3827 static const Real dzetadxi[4][3] =
3837 static const unsigned short int independent_var_indices[6][2] =
3851 static const unsigned short int zeta_indices[10][2] =
3866 const unsigned int my_j = independent_var_indices[j][0];
3867 const unsigned int my_k = independent_var_indices[j][1];
3871 return 4.*dzetadxi[i][my_j]*dzetadxi[i][my_k];
3876 const unsigned short int my_m = zeta_indices[i][0];
3877 const unsigned short int my_n = zeta_indices[i][1];
3879 return 4.*(dzetadxi[my_n][my_j]*dzetadxi[my_m][my_k] +
3880 dzetadxi[my_m][my_j]*dzetadxi[my_n][my_k] );
3883 libmesh_error_msg(
"Invalid shape function index " << i);
3891 libmesh_assert_less (i, 15);
3893 const Real xi = p(0);
3894 const Real eta = p(1);
3895 const Real zeta = p(2);
3906 return 2.*(1. - zeta);
3919 return 2.*(1. + zeta);
3921 return 4.*(zeta - 1);
3923 return -4.*(1. + zeta);
3925 libmesh_error_msg(
"Invalid i = " << i);
3936 return 2.*(1. - zeta);
3947 return 2.*(1. + zeta);
3950 return 2.*(zeta - 1.);
3953 return -2.*(1. + zeta);
3955 libmesh_error_msg(
"Invalid i = " << i);
3966 return 2.*(1. - zeta);
3979 return 2.*(1. + zeta);
3981 return 4.*(zeta - 1.);
3983 return -4.*(1. + zeta);
3985 libmesh_error_msg(
"Invalid i = " << i);
3995 return 1.5 - zeta - 2.*xi - 2.*eta;
3997 return 0.5 + zeta - 2.*xi;
4003 return -1.5 - zeta + 2.*xi + 2.*eta;
4005 return -0.5 + zeta + 2.*xi;
4007 return 4.*xi + 2.*eta - 2.;
4017 return -4.*xi - 2.*eta + 2.;
4023 libmesh_error_msg(
"Invalid i = " << i);
4033 return 1.5 - zeta - 2.*xi - 2.*eta;
4039 return .5 + zeta - 2.*eta;
4041 return -1.5 - zeta + 2.*xi + 2.*eta;
4043 return -.5 + zeta + 2.*eta;
4049 return 2.*xi + 4.*eta - 2.;
4059 return -2.*xi - 4.*eta + 2.;
4061 libmesh_error_msg(
"Invalid i = " << i);
4072 return 1. - xi - eta;
4087 return 2.*xi + 2.*eta - 2.;
4093 libmesh_error_msg(
"Invalid i = " << i);
4098 libmesh_error_msg(
"Invalid j = " << j);
4107 "High order on first order elements only supported for L2 families");
4108 libmesh_fallthrough();
4113 libmesh_assert_less (i, 18);
4118 Point p2d(p(0),p(1));
4122 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2};
4123 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2, 3, 4, 5, 0, 1, 2, 3, 4, 5, 3, 4, 5};
4158 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
4168 "High order on first order elements only supported for L2 families");
4169 libmesh_fallthrough();
4173 libmesh_assert_less (i, 14);
4175 const Real xi = p(0);
4176 const Real eta = p(1);
4177 const Real zeta = p(2);
4178 const Real eps = 1.e-35;
4183 p1 = 0.5*(1. - eta - zeta),
4184 p2 = 0.5*(1. + xi - zeta),
4185 p3 = 0.5*(1. + eta - zeta),
4186 p4 = 0.5*(1. - xi - zeta);
4191 den = (-1. + zeta + eps),
4198 numer_mp = xi*eta - zeta + zeta*zeta,
4199 numer_pm = xi*eta + zeta - zeta*zeta;
4209 return -p1*eta/den2;
4220 return 2.*p1*eta/den2;
4223 return 4.*p1*p3/den2;
4225 return -2.*p3*eta/den2;
4227 return -8.*p1*p3/den2;
4229 libmesh_error_msg(
"Invalid i = " << i);
4238 return 0.25*numer_mp/den2
4244 return 0.25*numer_pm/den2
4250 return 0.25*numer_mp/den2
4256 return 0.25*numer_pm/den2
4297 return 4.*p4*p1/den2
4303 libmesh_error_msg(
"Invalid i = " << i);
4326 return 4.*p2*p4/den2;
4328 return -2.*p2*xi/den2;
4330 return 2.*p4*xi/den2;
4332 return -8.*p2*p4/den2;
4334 libmesh_error_msg(
"Invalid i = " << i);
4344 return 0.25*numer_mp/den2
4345 - 0.5*p1*(2.*zeta - 1.)/den2
4349 - 2.*p4*p1*eta/den3;
4352 return 0.25*numer_pm/den2
4353 - 0.5*p1*(1 - 2.*zeta)/den2
4357 + 2.*p1*p2*eta/den3;
4360 return -0.25*numer_mp/den2
4361 + 0.5*p3*(2.*zeta - 1.)/den2
4365 - 2.*p2*p3*eta/den3;
4368 return -0.25*numer_pm/den2
4369 + 0.5*p3*(1 - 2.*zeta)/den2
4373 + 2.*p3*p4*eta/den3;
4382 - 4.*p1*p2*eta/den3;
4397 + 4.*p2*p3*eta/den3;
4429 return -4.*p4*p1/den2
4434 + 16.*p1*p2*p3/den3;
4437 libmesh_error_msg(
"Invalid i = " << i);
4446 return 0.25*numer_mp/den2
4447 - 0.5*p4*(2.*zeta - 1.)/den2
4454 return -0.25*numer_pm/den2
4455 + 0.5*p2*(1. - 2.*zeta)/den2
4462 return -0.25*numer_mp/den2
4463 + 0.5*p2*(2.*zeta - 1.)/den2
4470 return 0.25*numer_pm/den2
4471 - 0.5*p4*(1. - 2.*zeta)/den2
4531 return 4.*p3*p4/den2
4536 - 16.*p2*p1*p4/den3;
4539 libmesh_error_msg(
"Invalid i = " << i);
4548 return 0.5*numer_mp/den2
4549 - p1*(2.*zeta - 1.)/den2
4550 + 2.*p1*numer_mp/den3
4551 - p4*(2.*zeta - 1.)/den2
4552 + 2.*p4*numer_mp/den3
4554 - 4.*p4*p1*(2.*zeta - 1.)/den3
4555 + 6.*p4*p1*numer_mp/den4;
4558 return -0.5*numer_pm/den2
4559 + p2*(1 - 2.*zeta)/den2
4560 - 2.*p2*numer_pm/den3
4561 + p1*(1 - 2.*zeta)/den2
4562 - 2.*p1*numer_pm/den3
4564 + 4.*p1*p2*(1 - 2.*zeta)/den3
4565 - 6.*p1*p2*numer_pm/den4;
4568 return 0.5*numer_mp/den2
4569 - p3*(2.*zeta - 1.)/den2
4570 + 2.*p3*numer_mp/den3
4571 - p2*(2.*zeta - 1.)/den2
4572 + 2.*p2*numer_mp/den3
4574 - 4.*p2*p3*(2.*zeta - 1.)/den3
4575 + 6.*p2*p3*numer_mp/den4;
4578 return -0.5*numer_pm/den2
4579 + p4*(1 - 2.*zeta)/den2
4580 - 2.*p4*numer_pm/den3
4581 + p3*(1 - 2.*zeta)/den2
4582 - 2.*p3*numer_pm/den3
4584 + 4.*p3*p4*(1 - 2.*zeta)/den3
4585 - 6.*p3*p4*numer_pm/den4;
4591 return -2.*p1*eta/den2
4597 - 24.*p2*p1*p4*eta/den4;
4600 return 2.*p3*xi/den2
4606 + 24.*p1*p2*p3*xi/den4;
4609 return 2.*p4*eta/den2
4615 + 24.*p2*p3*p4*eta/den4;
4618 return -2.*p1*xi/den2
4624 - 24.*p3*p4*p1*xi/den4;
4633 - 8.*p1*p4*zeta/den3;
4642 - 8.*p2*p1*zeta/den3;
4651 - 8.*p3*p2*zeta/den3;
4660 - 8.*p4*p3*zeta/den3;
4663 return 8.*p3*p4/den2
4673 + 96.*p1*p2*p3*p4/den4;
4676 libmesh_error_msg(
"Invalid i = " << i);
4681 libmesh_error_msg(
"Invalid j = " << j);
4690 libmesh_assert_less (i, 13);
4692 const Real xi = p(0);
4693 const Real eta = p(1);
4694 const Real zeta = p(2);
4695 const Real eps = 1.e-35;
4700 den = (-1. + zeta + eps),
4716 return 0.5*(-1. + zeta + eta)/den;
4720 return 0.5*(-1. + zeta - eta)/den;
4732 return (1. - eta - zeta)/den;
4735 return (1. + eta - zeta)/den;
4738 libmesh_error_msg(
"Invalid i = " << i);
4747 return 0.25*(-1. + 2.*zeta + 2.*xi + 2.*eta)/den;
4750 return -0.25*(-1. + 2.*zeta - 2.*xi + 2.*eta)/den;
4753 return -0.25*(1. - 2.*zeta + 2.*xi + 2.*eta)/den;
4756 return 0.25*(1. - 2.*zeta - 2.*xi + 2.*eta)/den;
4786 libmesh_error_msg(
"Invalid i = " << i);
4797 return 0.5*(-1. + zeta + xi)/den;
4801 return 0.5*(-1. + zeta - xi)/den;
4813 return (1. + xi - zeta)/den;
4816 return (1. - xi - zeta)/den;
4819 libmesh_error_msg(
"Invalid i = " << i);
4829 return -0.25*(-1. + 2.*zeta - zeta2 + eta + 2.*eta*xi + eta2)/den2;
4832 return 0.25*(-1. + 2.*zeta - zeta2 + eta - 2.*eta*xi + eta2)/den2;
4835 return 0.25*(-1. + 2.*zeta - zeta2 - eta + 2.*eta*xi + eta2)/den2;
4838 return -0.25*(-1. + 2.*zeta - zeta2 - eta - 2.*eta*xi + eta2)/den2;
4847 return -0.5*(1. + zeta2 + eta2 - 2.*zeta)/den2;
4850 return -eta*xi/den2;
4853 return 0.5*(1. + zeta2 + eta2 - 2.*zeta)/den2;
4856 return (-1. - zeta2 + eta + 2.*zeta)/den2;
4859 return -(-1. - zeta2 + eta + 2.*zeta)/den2;
4862 return (1. + zeta2 + eta - 2.*zeta)/den2;
4865 return -(1. + zeta2 + eta - 2.*zeta)/den2;
4868 libmesh_error_msg(
"Invalid i = " << i);
4877 return -0.25*(-1. + 2.*zeta - zeta2 + xi + 2.*eta*xi + xi2)/den2;
4880 return 0.25*(1. - 2.*zeta + zeta2 + xi + 2.*eta*xi - xi2)/den2;
4883 return 0.25*(-1. + 2.*zeta - zeta2 - xi + 2.*eta*xi + xi2)/den2;
4886 return -0.25*(1. - 2.*zeta + zeta2 - xi + 2.*eta*xi - xi2)/den2;
4892 return 0.5*(1. + xi2 + zeta2 - 2.*zeta)/den2;
4895 return -eta*xi/den2;
4898 return -0.5*(1. + xi2 + zeta2 - 2.*zeta)/den2;
4904 return (-1. - zeta2 + xi + 2.*zeta)/den2;
4907 return -(1. + zeta2 + xi - 2.*zeta)/den2;
4910 return (1. + zeta2 + xi - 2.*zeta)/den2;
4913 return -(-1. - zeta2 + xi + 2.*zeta)/den2;
4916 libmesh_error_msg(
"Invalid i = " << i);
4925 return 0.5*(xi + eta + 1.)*eta*xi/den3;
4928 return -0.5*(eta - xi + 1.)*eta*xi/den3;
4931 return -0.5*(xi + eta - 1.)*eta*xi/den3;
4934 return 0.5*(eta - xi - 1.)*eta*xi/den3;
4940 return -(1. - 3.*zeta + 3.*zeta2 - zeta3 + eta*xi2)/den3;
4943 return (-1. + 3.*zeta - 3.*zeta2 + zeta3 + eta2*xi)/den3;
4946 return (-1. + 3.*zeta - 3.*zeta2 + zeta3 + eta*xi2)/den3;
4949 return -(1. - 3.*zeta + 3.*zeta2 - zeta3 + eta2*xi)/den3;
4952 return -2.*(-1. + 3.*zeta - 3.*zeta2 + zeta3 + eta*xi)/den3;
4955 return 2.*(1. - 3.*zeta + 3.*zeta2 - zeta3 + eta*xi)/den3;
4958 return -2.*(-1. + 3.*zeta - 3.*zeta2 + zeta3 + eta*xi)/den3;
4961 return 2.*(1. - 3.*zeta + 3.*zeta2 - zeta3 + eta*xi)/den3;
4964 libmesh_error_msg(
"Invalid i = " << i);
4969 libmesh_error_msg(
"Invalid j = " << j);
4985 libmesh_assert_less (i, 14);
4993 static const Real dzetadxi[4][3] =
5003 static const unsigned short int independent_var_indices[6][2] =
5017 static const unsigned short int zeta_indices[10][2] =
5032 const unsigned int my_j = independent_var_indices[j][0];
5033 const unsigned int my_k = independent_var_indices[j][1];
5037 returnval = 4.*dzetadxi[i][my_j]*dzetadxi[i][my_k];
5041 const unsigned short int my_m = zeta_indices[i][0];
5042 const unsigned short int my_n = zeta_indices[i][1];
5045 4.*(dzetadxi[my_n][my_j]*dzetadxi[my_m][my_k] +
5046 dzetadxi[my_m][my_j]*dzetadxi[my_n][my_k] );
5049 const Real zeta1 = p(0);
5050 const Real zeta2 = p(1);
5051 const Real zeta3 = p(2);
5052 const Real zeta0 = 1. - zeta1 - zeta2 - zeta3;
5056 Real d2bubble012, d2bubble013, d2bubble023, d2bubble123;
5062 d2bubble012 = -2.*zeta2;
5063 d2bubble013 = -2.*zeta3;
5072 d2bubble012 = (zeta0-zeta1)-zeta2;
5073 d2bubble013 = -zeta3;
5074 d2bubble123 = zeta3;
5075 d2bubble023 = -zeta3;
5082 d2bubble012 = -2.*zeta1;
5085 d2bubble023 = -2.*zeta3;
5092 d2bubble012 = -zeta2;
5093 d2bubble013 = (zeta0-zeta3)-zeta1;
5094 d2bubble123 = zeta2;
5095 d2bubble023 = -zeta2;
5102 d2bubble012 = -zeta1;
5103 d2bubble013 = -zeta1;
5104 d2bubble123 = zeta1;
5105 d2bubble023 = (zeta0-zeta3)-zeta2;
5113 d2bubble013 = -2.*zeta1;
5115 d2bubble023 = -2.*zeta2;
5120 libmesh_error_msg(
"Invalid j = " << j);
5126 return returnval + 3.*(d2bubble012+d2bubble013+d2bubble023);
5129 return returnval + 3.*(d2bubble012+d2bubble013+d2bubble123);
5132 return returnval + 3.*(d2bubble012+d2bubble023+d2bubble123);
5135 return returnval + 3.*(d2bubble013+d2bubble023+d2bubble123);
5138 return returnval - 12.*(d2bubble012+d2bubble013);
5141 return returnval - 12.*(d2bubble012+d2bubble123);
5144 return returnval - 12.*(d2bubble012+d2bubble023);
5147 return returnval - 12.*(d2bubble013+d2bubble023);
5150 return returnval - 12.*(d2bubble013+d2bubble123);
5153 return returnval - 12.*(d2bubble023+d2bubble123);
5156 return 27.*d2bubble012;
5159 return 27.*d2bubble013;
5162 return 27.*d2bubble123;
5165 return 27.*d2bubble023;
5168 libmesh_error_msg(
"Invalid i = " << i);
5174 libmesh_assert_less (i, 20);
5179 Point p2d(p(0),p(1));
5183 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2, 0, 1};
5237 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
5241 return mainval - bubbleval / 9;
5243 return mainval + bubbleval * (
Real(4) / 9);
5248 libmesh_assert_less (i, 21);
5253 Point p2d(p(0),p(1));
5257 static const unsigned int i0[] = {0, 0, 0, 1, 1, 1, 0, 0, 0, 2, 2, 2, 1, 1, 1, 2, 2, 2, 0, 1, 2};
5258 static const unsigned int i1[] = {0, 1, 2, 0, 1, 2, 3, 4, 5, 0, 1, 2, 3, 4, 5, 3, 4, 5, 6, 6, 6};
5293 libmesh_error_msg(
"Invalid shape function derivative j = " << j);
5300 fe_lagrange_3D_shape_deriv<T>);
5310 libmesh_error_msg(
"ERROR: Unsupported 3D FE order!: " << order);
5313 #else // LIBMESH_DIM != 3 5315 libmesh_not_implemented();
5319 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES Real fe_lagrange_1D_quadratic_shape_second_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real)
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Real fe_lagrange_1D_linear_shape_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real)
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Real fe_lagrange_1D_quadratic_shape_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real xi)
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
Real fe_lagrange_1D_quadratic_shape(const unsigned int i, const Real xi)
static void all_shape_derivs(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> *comps[3], const bool add_p_level=true)
Fills comps with dphidxi (and in higher dimensions, eta/zeta) derivative component values for all sha...
static void default_all_shape_derivs(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> *comps[3], const bool add_p_level=true)
A default implementation for all_shape_derivs.
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.
static void default_all_shapes(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> &v, const bool add_p_level=true)
A default implementation for all_shapes.
OutputShape fe_fdm_deriv(const Elem *elem, const Order order, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level, OutputShape(*shape_func)(const Elem *, const Order, const unsigned int, const Point &, const bool))
Helper functions for finite differenced derivatives in cases where analytical calculations haven't be...
LIBMESH_DEFAULT_VECTORIZED_FE(template<>Real FE< 0, BERNSTEIN)
A specific instantiation of the FEBase class.
void libmesh_ignore(const Args &...)
std::tuple< unsigned int, Real, Real, Real > subelement_coordinates(const Point &p, Real tol=TOLERANCE *TOLERANCE) const
OutputShape fe_fdm_second_deriv(const Elem *elem, const Order order, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level, OutputShape(*deriv_func)(const Elem *, const Order, const unsigned int, const unsigned int, const Point &, const bool))
static void default_shape_derivs(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
A default implementation for shape_derivs.
std::string enum_to_string(const T e)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static void shape_derivs(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
Fills v with the derivative of the shape function, evaluated at all points p.
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
static void shapes(const Elem *elem, const Order o, const unsigned int i, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
Fills v with the values of the shape function, evaluated at all points p.
The C0Polyhedron is an element in 3D with an arbitrary (but fixed) number of polygonal first-order (C...
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
This class defines a tensor in LIBMESH_DIM dimensional Real or Complex space.
static void default_shapes(const Elem *elem, const Order o, const unsigned int i, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
A default implementation for shapes.
Real fe_lagrange_1D_linear_shape(const unsigned int i, const Real xi)
static void all_shapes(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape > > &v, const bool add_p_level=true)
Fills v[i][qp] with the values of the shape functions, evaluated at all points in p...