20 #include "libmesh/fe.h" 21 #include "libmesh/elem.h" 22 #include "libmesh/number_lookups.h" 23 #include "libmesh/enum_to_string.h" 32 Real fe_triangle_helper (
const Elem & elem,
33 const Real edgenumerator,
35 const unsigned int basisorder,
36 const Order totalorder,
37 const unsigned int noden);
40 Real fe_hierarchic_2D_shape(
const Elem * elem,
44 const bool add_p_level);
47 Real fe_hierarchic_2D_shape_deriv(
const Elem * elem,
52 const bool add_p_level);
54 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 57 Real fe_hierarchic_2D_shape_second_deriv(
const Elem * elem,
62 const bool add_p_level);
64 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 67 std::tuple<unsigned int, unsigned int, Real>
68 quad_indices(
const Elem * elem,
69 const unsigned int totalorder,
72 libmesh_assert_less (i, (totalorder+1u)*(totalorder+1u));
91 else if (i < totalorder + 3u)
92 { i0 = i - 2; i1 = 0; }
93 else if (i < 2u*totalorder + 2)
94 { i0 = 1; i1 = i - totalorder - 1; }
95 else if (i < 3u*totalorder + 1)
96 { i0 = i - 2u*totalorder; i1 = 1; }
97 else if (i < 4u*totalorder)
98 { i0 = 0; i1 = i - 3u*totalorder + 1; }
102 unsigned int basisnum = i - 4*totalorder;
111 if ((i0%2) && (i0 > 2) && (i1 == 0))
113 else if ((i0%2) && (i0>2) && (i1 == 1))
115 else if ((i0 == 0) && (i1%2) && (i1>2))
117 else if ((i0 == 1) && (i1%2) && (i1>2))
142 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
154 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
166 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
175 const unsigned int i,
177 const bool add_p_level)
179 return fe_hierarchic_2D_shape<HIERARCHIC>(elem, order, i, p, add_p_level);
187 const unsigned int i,
189 const bool add_p_level)
191 return fe_hierarchic_2D_shape<HIERARCHIC>(elem, fet.
order, i, p, add_p_level);
198 const unsigned int i,
200 const bool add_p_level)
202 return fe_hierarchic_2D_shape<L2_HIERARCHIC>(elem, order, i, p, add_p_level);
209 const unsigned int i,
211 const bool add_p_level)
213 return fe_hierarchic_2D_shape<L2_HIERARCHIC>(elem, fet.
order, i, p, add_p_level);
220 const unsigned int i,
222 const bool add_p_level)
227 const Order totalorder = order + add_p_level*elem->
p_level();
229 const unsigned int dofs_per_side = totalorder+1u;
236 libmesh_assert_less(i, 3*dofs_per_side);
244 const Real zeta1 = p(0);
245 const Real zeta2 = p(1);
246 const Real zeta0 = 1. - zeta1 - zeta2;
248 if (zeta1 > zeta2 && zeta0 > zeta2)
250 if (i >= dofs_per_side)
256 if ((i < 2 || i % 2) &&
262 else if (zeta1 > zeta0 && zeta2 > zeta0)
264 if (i < dofs_per_side ||
265 i >= 2*dofs_per_side)
271 const unsigned int side_i = i - dofs_per_side;
273 if ((side_i < 2 || side_i % 2) &&
283 if (i < 2*dofs_per_side)
289 const unsigned int side_i = i - 2*dofs_per_side;
291 if ((side_i < 2 || side_i % 2) &&
303 libmesh_assert_less(i, 4*dofs_per_side);
311 const Real xi = p(0), eta = p(1);
316 if (i >= dofs_per_side)
322 if ((i < 2 || i % 2) &&
330 if (i < dofs_per_side ||
331 i >= 2*dofs_per_side)
337 const unsigned int side_i = i - dofs_per_side;
339 if ((side_i < 2 || side_i % 2) &&
350 if (i < 2*dofs_per_side ||
351 i >= 3*dofs_per_side)
357 const unsigned int side_i = i - 2*dofs_per_side;
359 if ((side_i < 2 || side_i % 2) &&
367 if (i < 3*dofs_per_side)
373 const unsigned int side_i = i - 3*dofs_per_side;
375 if ((side_i < 2 || side_i % 2) &&
393 const unsigned int i,
395 const bool add_p_level)
408 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
421 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
434 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
443 const unsigned int i,
444 const unsigned int j,
446 const bool add_p_level)
448 return fe_hierarchic_2D_shape_deriv<HIERARCHIC>(elem, order, i, j, p, add_p_level);
455 const unsigned int i,
456 const unsigned int j,
458 const bool add_p_level)
460 return fe_hierarchic_2D_shape_deriv<HIERARCHIC>(elem, fet.
order, i, j, p, add_p_level);
469 const unsigned int i,
470 const unsigned int j,
472 const bool add_p_level)
474 return fe_hierarchic_2D_shape_deriv<L2_HIERARCHIC>(elem, order, i, j, p, add_p_level);
481 const unsigned int i,
482 const unsigned int j,
484 const bool add_p_level)
486 return fe_hierarchic_2D_shape_deriv<L2_HIERARCHIC>(elem, fet.
order, i, j, p, add_p_level);
494 const unsigned int i,
495 const unsigned int j,
497 const bool add_p_level)
503 const Order totalorder = order + add_p_level*elem->
p_level();
508 const unsigned int dofs_per_side = totalorder+1u;
519 libmesh_assert_less(i, 3*dofs_per_side);
520 libmesh_assert_less (j, 2);
528 const Real zeta1 = p(0);
529 const Real zeta2 = p(1);
530 const Real zeta0 = 1. - zeta1 - zeta2;
532 if (zeta1 > zeta2 && zeta0 > zeta2)
537 if (i >= dofs_per_side)
543 if ((i < 2 || i % 2) &&
549 else if (zeta1 > zeta0 && zeta2 > zeta0)
551 if (i < dofs_per_side ||
552 i >= 2*dofs_per_side)
558 const unsigned int side_i = i - dofs_per_side;
560 if ((side_i < 2 || side_i % 2) &&
577 if (i < 2*dofs_per_side)
583 const unsigned int side_i = i - 2*dofs_per_side;
585 if ((side_i < 2 || side_i % 2) &&
598 libmesh_assert_less(i, 4*dofs_per_side);
606 const Real xi = p(0), eta = p(1);
611 if (i >= dofs_per_side)
615 if ((i < 2 || i % 2) &&
623 if (i < dofs_per_side ||
624 i >= 2*dofs_per_side)
629 const unsigned int side_i = i - dofs_per_side;
631 if ((side_i < 2 || side_i % 2) &&
642 if (i < 2*dofs_per_side ||
643 i >= 3*dofs_per_side)
648 const unsigned int side_i = i - 2*dofs_per_side;
650 if ((side_i < 2 || side_i % 2) &&
658 if (i < 3*dofs_per_side)
663 const unsigned int side_i = i - 3*dofs_per_side;
665 if ((side_i < 2 || side_i % 2) &&
683 const unsigned int i,
684 const unsigned int j,
686 const bool add_p_level)
692 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 701 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
714 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
727 libmesh_error_msg(
"Hierarchic shape functions require an Elem for edge orientation.");
736 const unsigned int i,
737 const unsigned int j,
739 const bool add_p_level)
741 return fe_hierarchic_2D_shape_second_deriv<HIERARCHIC>(elem, order, i, j, p, add_p_level);
748 const unsigned int i,
749 const unsigned int j,
751 const bool add_p_level)
753 return fe_hierarchic_2D_shape_second_deriv<HIERARCHIC>(elem, fet.
order, i, j, p, add_p_level);
760 const unsigned int i,
761 const unsigned int j,
763 const bool add_p_level)
765 return fe_hierarchic_2D_shape_second_deriv<L2_HIERARCHIC>(elem, order, i, j, p, add_p_level);
772 const unsigned int i,
773 const unsigned int j,
775 const bool add_p_level)
777 return fe_hierarchic_2D_shape_second_deriv<L2_HIERARCHIC>(elem, fet.
order, i, j, p, add_p_level);
784 const unsigned int i,
785 const unsigned int j,
787 const bool add_p_level)
792 const Order totalorder = order + add_p_level*elem->
p_level();
797 const unsigned int dofs_per_side = totalorder+1u;
812 libmesh_assert_less(i, 4*dofs_per_side);
820 const Real xi = p(0), eta = p(1);
825 if (i >= dofs_per_side)
829 if ((i < 2 || i % 2) &&
837 if (i < dofs_per_side ||
838 i >= 2*dofs_per_side)
843 const unsigned int side_i = i - dofs_per_side;
845 if ((side_i < 2 || side_i % 2) &&
856 if (i < 2*dofs_per_side ||
857 i >= 3*dofs_per_side)
862 const unsigned int side_i = i - 2*dofs_per_side;
864 if ((side_i < 2 || side_i % 2) &&
872 if (i < 3*dofs_per_side)
877 const unsigned int side_i = i - 3*dofs_per_side;
879 if ((side_i < 2 || side_i % 2) &&
897 const unsigned int i,
898 const unsigned int j,
900 const bool add_p_level)
905 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES 915 Real fe_triangle_helper (
const Elem & elem,
916 const Real edgenumerator,
918 const unsigned int basisorder,
919 const Order totalorder,
920 const unsigned int noden)
924 if (basisorder%2 && (elem.
point(noden) > elem.
point((noden+1)%3)))
933 unsigned int basisfactorial = 1.;
934 for (
unsigned int n=2; n <= basisorder; ++n)
937 return std::pow(edgenumerator, basisorder) / basisfactorial;
942 const Real edgeval = edgenumerator / crossval;
945 return flip * crossfunc *
947 basisorder, edgeval);
950 template <FEFamily T>
951 Real fe_hierarchic_2D_shape(
const Elem * elem,
953 const unsigned int i,
955 const bool add_p_level)
959 const Order totalorder = order + add_p_level*elem->
p_level();
960 libmesh_assert_greater (totalorder, 0);
962 switch (elem->
type())
969 const Real zeta1 = p(0);
970 const Real zeta2 = p(1);
971 const Real zeta0 = 1. - zeta1 - zeta2;
973 libmesh_assert_less (i, (totalorder+1u)*(totalorder+2u)/2);
975 elem->
type() ==
TRI7 || totalorder < 2);
985 else if (i < totalorder + 2u)
987 const unsigned int basisorder = i - 1;
989 const Real crossval = zeta0 + zeta1;
990 const Real edgenumerator = zeta1 - zeta0;
992 return fe_triangle_helper(*elem, edgenumerator, crossval,
993 basisorder, totalorder, 0);
995 else if (i < 2u*totalorder + 1)
997 const unsigned int basisorder = i - totalorder;
999 const Real crossval = zeta2 + zeta1;
1000 const Real edgenumerator = zeta2 - zeta1;
1002 return fe_triangle_helper(*elem, edgenumerator, crossval,
1003 basisorder, totalorder, 1);
1005 else if (i < 3u*totalorder)
1007 const unsigned int basisorder = i - (2u*totalorder) + 1;
1009 const Real crossval = zeta0 + zeta2;
1010 const Real edgenumerator = zeta0 - zeta2;
1012 return fe_triangle_helper(*elem, edgenumerator, crossval,
1013 basisorder, totalorder, 2);
1018 const unsigned int basisnum = i - (3u*totalorder);
1024 for (
unsigned int n = 0; n != exp0; ++n)
1026 for (
unsigned int n = 0; n != exp1; ++n)
1037 libmesh_fallthrough();
1044 auto [i0, i1, f] = quad_indices(elem, totalorder, i);
1058 template <FEFamily T>
1059 Real fe_hierarchic_2D_shape_deriv(
const Elem * elem,
1061 const unsigned int i,
1062 const unsigned int j,
1064 const bool add_p_level)
1070 const Order totalorder = order + add_p_level*elem->
p_level();
1072 libmesh_assert_greater (totalorder, 0);
1088 libmesh_fallthrough();
1095 auto [i0, i1, f] = quad_indices(elem, totalorder, i);
1110 libmesh_error_msg(
"Invalid derivative index j = " << j);
1123 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 1125 template <FEFamily T>
1126 Real fe_hierarchic_2D_shape_second_deriv(
const Elem * elem,
1128 const unsigned int i,
1129 const unsigned int j,
1131 const bool add_p_level)
1137 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
This is the base class from which all geometric element types are derived.
const unsigned char triangular_number_row[]
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.
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...
const unsigned char square_number_column[]
LIBMESH_DEFAULT_VECTORIZED_FE(template<>Real FE< 0, BERNSTEIN)
A specific instantiation of the FEBase class.
bool positive_edge_orientation(const unsigned int i) 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))
std::string enum_to_string(const T e)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const unsigned char triangular_number_column[]
const unsigned char square_number_row[]
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
const Point & point(const unsigned int i) const
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)