22 #include "libmesh/fe.h" 
   23 #include "libmesh/elem.h" 
   24 #include "libmesh/fe_interface.h" 
   37 static const Elem * old_elem_ptr = 
nullptr;
 
   42 static Real d1xd1x, d2xd2x;
 
   44 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
   46 Real clough_raw_shape_second_deriv(
const unsigned int basis_num,
 
   47                                    const unsigned int deriv_type,
 
   52 Real clough_raw_shape_deriv(
const unsigned int basis_num,
 
   53                             const unsigned int deriv_type,
 
   55 Real clough_raw_shape(
const unsigned int basis_num,
 
   60 void clough_compute_coefs(
const Elem * elem)
 
   69   if (elem->
id() == old_elem_id &&
 
   76   const Order mapping_order        (elem->default_order());
 
   77   const ElemType mapping_elem_type (elem->type());
 
   79   const FEType map_fe_type(mapping_order, mapping_family);
 
   81   const int n_mapping_shape_functions =
 
   85   std::vector<Point> dofpt;
 
   86   dofpt.push_back(
Point(0));
 
   87   dofpt.push_back(
Point(1));
 
   90   std::vector<Real> dxdxi(2);
 
   91   std::vector<Real> dxidx(2);
 
   96   for (
int p = 0; p != 2; ++p)
 
   98       for (
int i = 0; i != n_mapping_shape_functions; ++i)
 
  100           const Real ddxi = shape_deriv_ptr
 
  101             (elem, mapping_order, i, 0, dofpt[p], 
false);
 
  102           dxdxi[p] += dofpt[p](0) * ddxi;
 
  110   if (elem->id() == old_elem_id &&
 
  111       elem == old_elem_ptr)
 
  113       libmesh_assert_equal_to(d1xd1x, dxdxi[0]);
 
  114       libmesh_assert_equal_to(d2xd2x, dxdxi[1]);
 
  118   old_elem_id = elem->id();
 
  126 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  129 Real clough_raw_shape_second_deriv(
const unsigned int basis_num,
 
  130                                    const unsigned int deriv_type,
 
  152           libmesh_error_msg(
"Invalid shape function index i = " <<
 
  157       libmesh_error_msg(
"Invalid shape function derivative j = " <<
 
  165 Real clough_raw_shape_deriv(
const unsigned int basis_num,
 
  166                             const unsigned int deriv_type,
 
  177           return -6*xi + 6*xi*xi;
 
  179           return 6*xi - 6*xi*xi;
 
  181           return 1 - 4*xi + 3*xi*xi;
 
  183           return -2*xi + 3*xi*xi;
 
  186           libmesh_error_msg(
"Invalid shape function index i = " <<
 
  191       libmesh_error_msg(
"Invalid shape function derivative j = " <<
 
  196 Real clough_raw_shape(
const unsigned int basis_num,
 
  204       return 1 - 3*xi*xi + 2*xi*xi*xi;
 
  206       return 3*xi*xi - 2*xi*xi*xi;
 
  208       return xi - 2*xi*xi + xi*xi*xi;
 
  210       return -xi*xi + xi*xi*xi;
 
  213       libmesh_error_msg(
"Invalid shape function index i = " <<
 
  232   libmesh_error_msg(
"Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
 
  241                          const unsigned int i,
 
  243                          const bool add_p_level)
 
  247   clough_compute_coefs(elem);
 
  251   const Order totalorder =
 
  252     static_cast<Order>(order + add_p_level * elem->
p_level());
 
  265               libmesh_assert_less (i, 4);
 
  270                   return clough_raw_shape(0, p);
 
  272                   return clough_raw_shape(1, p);
 
  274                   return d1xd1x * clough_raw_shape(2, p);
 
  276                   return d2xd2x * clough_raw_shape(3, p);
 
  278                   libmesh_error_msg(
"Invalid shape function index i = " << i);
 
  282             libmesh_error_msg(
"ERROR: Unsupported element type = " << type);
 
  287       libmesh_error_msg(
"ERROR: Unsupported polynomial order = " << totalorder);
 
  300   libmesh_error_msg(
"Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
 
  309                                const unsigned int i,
 
  310                                const unsigned int j,
 
  312                                const bool add_p_level)
 
  316   clough_compute_coefs(elem);
 
  320   const Order totalorder =
 
  321     static_cast<Order>(order + add_p_level * elem->
p_level());
 
  337                   return clough_raw_shape_deriv(0, j, p);
 
  339                   return clough_raw_shape_deriv(1, j, p);
 
  341                   return d1xd1x * clough_raw_shape_deriv(2, j, p);
 
  343                   return d2xd2x * clough_raw_shape_deriv(3, j, p);
 
  345                   libmesh_error_msg(
"Invalid shape function index i = " << i);
 
  349             libmesh_error_msg(
"ERROR: Unsupported element type = " << type);
 
  354       libmesh_error_msg(
"ERROR: Unsupported polynomial order = " << totalorder);
 
  359 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  368   libmesh_error_msg(
"Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
 
  377                                       const unsigned int i,
 
  378                                       const unsigned int j,
 
  380                                       const bool add_p_level)
 
  384   clough_compute_coefs(elem);
 
  388   const Order totalorder =
 
  389     static_cast<Order>(order + add_p_level * elem->
p_level());
 
  405                   return clough_raw_shape_second_deriv(0, j, p);
 
  407                   return clough_raw_shape_second_deriv(1, j, p);
 
  409                   return d1xd1x * clough_raw_shape_second_deriv(2, j, p);
 
  411                   return d2xd2x * clough_raw_shape_second_deriv(3, j, p);
 
  413                   libmesh_error_msg(
"Invalid shape function index i = " << i);
 
  417             libmesh_error_msg(
"ERROR: Unsupported element type = " << type);
 
  422       libmesh_error_msg(
"ERROR: Unsupported polynomial order = " << totalorder);