36 static const Elem * old_elem_ptr =
nullptr;
41 static Real d1xd1x, d2xd2x;
43 Real clough_raw_shape_second_deriv(
const unsigned int basis_num,
44 const unsigned int deriv_type,
46 Real clough_raw_shape_deriv(
const unsigned int basis_num,
47 const unsigned int deriv_type,
49 Real clough_raw_shape(
const unsigned int basis_num,
54 void clough_compute_coefs(
const Elem * elem)
63 if (elem->
id() == old_elem_id &&
68 const Order mapping_order (elem->default_order());
69 const ElemType mapping_elem_type (elem->type());
70 const int n_mapping_shape_functions =
75 std::vector<Point> dofpt;
76 dofpt.push_back(
Point(0));
77 dofpt.push_back(
Point(1));
80 std::vector<Real> dxdxi(2);
81 std::vector<Real> dxidx(2);
83 for (
int p = 0; p != 2; ++p)
85 for (
int i = 0; i != n_mapping_shape_functions; ++i)
88 (mapping_elem_type, mapping_order, i, 0, dofpt[p]);
89 dxdxi[p] += dofpt[p](0) * ddxi;
97 if (elem->id() == old_elem_id &&
100 libmesh_assert_equal_to(d1xd1x, dxdxi[0]);
101 libmesh_assert_equal_to(d2xd2x, dxdxi[1]);
105 old_elem_id = elem->id();
114 Real clough_raw_shape_second_deriv(
const unsigned int basis_num,
115 const unsigned int deriv_type,
137 libmesh_error_msg(
"Invalid shape function index i = " <<
142 libmesh_error_msg(
"Invalid shape function derivative j = " <<
149 Real clough_raw_shape_deriv(
const unsigned int basis_num,
150 const unsigned int deriv_type,
161 return -6*xi + 6*xi*xi;
163 return 6*xi - 6*xi*xi;
165 return 1 - 4*xi + 3*xi*xi;
167 return -2*xi + 3*xi*xi;
170 libmesh_error_msg(
"Invalid shape function index i = " <<
175 libmesh_error_msg(
"Invalid shape function derivative j = " <<
180 Real clough_raw_shape(
const unsigned int basis_num,
188 return 1 - 3*xi*xi + 2*xi*xi*xi;
190 return 3*xi*xi - 2*xi*xi*xi;
192 return xi - 2*xi*xi + xi*xi*xi;
194 return -xi*xi + xi*xi*xi;
197 libmesh_error_msg(
"Invalid shape function index i = " <<
216 libmesh_error_msg(
"Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
225 const unsigned int i,
228 libmesh_assert(elem);
230 clough_compute_coefs(elem);
247 libmesh_assert_less (i, 4);
252 return clough_raw_shape(0, p);
254 return clough_raw_shape(1, p);
256 return d1xd1x * clough_raw_shape(2, p);
258 return d2xd2x * clough_raw_shape(3, p);
260 libmesh_error_msg(
"Invalid shape function index i = " << i);
264 libmesh_error_msg(
"ERROR: Unsupported element type = " << type);
269 libmesh_error_msg(
"ERROR: Unsupported polynomial order = " << totalorder);
282 libmesh_error_msg(
"Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
291 const unsigned int i,
292 const unsigned int j,
295 libmesh_assert(elem);
297 clough_compute_coefs(elem);
317 return clough_raw_shape_deriv(0, j, p);
319 return clough_raw_shape_deriv(1, j, p);
321 return d1xd1x * clough_raw_shape_deriv(2, j, p);
323 return d2xd2x * clough_raw_shape_deriv(3, j, p);
325 libmesh_error_msg(
"Invalid shape function index i = " << i);
329 libmesh_error_msg(
"ERROR: Unsupported element type = " << type);
334 libmesh_error_msg(
"ERROR: Unsupported polynomial order = " << totalorder);
343 const unsigned int i,
344 const unsigned int j,
347 libmesh_assert(elem);
349 clough_compute_coefs(elem);
369 return clough_raw_shape_second_deriv(0, j, p);
371 return clough_raw_shape_second_deriv(1, j, p);
373 return d1xd1x * clough_raw_shape_second_deriv(2, j, p);
375 return d2xd2x * clough_raw_shape_second_deriv(3, j, p);
377 libmesh_error_msg(
"Invalid shape function index i = " << i);
381 libmesh_error_msg(
"ERROR: Unsupported element type = " << type);
386 libmesh_error_msg(
"ERROR: Unsupported polynomial order = " << totalorder);
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
virtual unsigned int n_shape_functions() const override
static const dof_id_type invalid_id
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)