libMesh
fe_clough_shape_1D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // libmesh includes
20 #include "libmesh/fe.h"
21 #include "libmesh/elem.h"
22 #include "libmesh/fe_interface.h"
23 #include "libmesh/enum_to_string.h"
24 
25 // Anonymous namespace for persistent variables.
26 // This allows us to cache the global-to-local mapping transformation
27 // This should also screw up multithreading royally
28 namespace
29 {
30 using namespace libMesh;
31 
32 // Keep track of which element was most recently used to generate
33 // cached data
34 static dof_id_type old_elem_id = DofObject::invalid_id;
35 static const Elem * old_elem_ptr = nullptr;
36 
37 // Coefficient naming: d(1)d(2n) is the coefficient of the
38 // global shape function corresponding to value 1 in terms of the
39 // local shape function corresponding to normal derivative 2
40 static Real d1xd1x, d2xd2x;
41 
42 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
43 
44 Real clough_raw_shape_second_deriv(const unsigned int basis_num,
45  const unsigned int deriv_type,
46  const Point & p);
47 
48 #endif
49 
50 Real clough_raw_shape_deriv(const unsigned int basis_num,
51  const unsigned int deriv_type,
52  const Point & p);
53 Real clough_raw_shape(const unsigned int basis_num,
54  const Point & p);
55 
56 
57 // Compute the static coefficients for an element
58 void clough_compute_coefs(const Elem * elem)
59 {
60  // Using static globals for old_elem_id, etc. will fail
61  // horribly with more than one thread.
62  libmesh_assert_equal_to (libMesh::n_threads(), 1);
63 
64  // Coefficients are cached from old elements; we rely on that cache
65  // except in dbg mode
66 #ifndef DEBUG
67  if (elem->id() == old_elem_id &&
68  elem == old_elem_ptr)
69  return;
70 #endif
71 
72  const FEFamily mapping_family = FEMap::map_fe_type(*elem);
73  const FEType map_fe_type(elem->default_order(), mapping_family);
74 
75  // Note: we explicitly don't consider the elem->p_level() when
76  // computing the number of mapping shape functions.
77  const int n_mapping_shape_functions =
78  FEInterface::n_shape_functions(map_fe_type, /*extra_order=*/0, elem);
79 
80  // Degrees of freedom are at vertices and edge midpoints
81  std::vector<Point> dofpt;
82  dofpt.push_back(Point(0));
83  dofpt.push_back(Point(1));
84 
85  // Mapping functions - first derivatives at each dofpt
86  std::vector<Real> dxdxi(2);
87  std::vector<Real> dxidx(2);
88 
89  FEInterface::shape_deriv_ptr shape_deriv_ptr =
90  FEInterface::shape_deriv_function(map_fe_type, elem);
91 
92  for (int p = 0; p != 2; ++p)
93  {
94  for (int i = 0; i != n_mapping_shape_functions; ++i)
95  {
96  const Real ddxi = shape_deriv_ptr
97  (map_fe_type, elem, i, 0, dofpt[p], /*add_p_level=*/false);
98  dxdxi[p] += dofpt[p](0) * ddxi;
99  }
100  }
101 
102  // Calculate derivative scaling factors
103 
104 #ifdef DEBUG
105  // The cached factors should equal our calculations
106  if (elem->id() == old_elem_id &&
107  elem == old_elem_ptr)
108  {
109  libmesh_assert_equal_to(d1xd1x, dxdxi[0]);
110  libmesh_assert_equal_to(d2xd2x, dxdxi[1]);
111  }
112 #endif
113 
114  old_elem_id = elem->id();
115  old_elem_ptr = elem;
116 
117  d1xd1x = dxdxi[0];
118  d2xd2x = dxdxi[1];
119 }
120 
121 
122 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
123 
124 // Return shape function second derivatives on the unit interval
125 Real clough_raw_shape_second_deriv(const unsigned int basis_num,
126  const unsigned int deriv_type,
127  const Point & p)
128 {
129  Real xi = p(0);
130 
131  switch (deriv_type)
132  {
133 
134  // second derivative in xi-xi direction
135  case 0:
136  switch (basis_num)
137  {
138  case 0:
139  return -6 + 12*xi;
140  case 1:
141  return 6 - 12*xi;
142  case 2:
143  return -4 + 6*xi;
144  case 3:
145  return -2 + 6*xi;
146 
147  default:
148  libmesh_error_msg("Invalid shape function index i = " <<
149  basis_num);
150  }
151 
152  default:
153  libmesh_error_msg("Invalid shape function derivative j = " <<
154  deriv_type);
155  }
156 }
157 
158 #endif
159 
160 
161 Real clough_raw_shape_deriv(const unsigned int basis_num,
162  const unsigned int deriv_type,
163  const Point & p)
164 {
165  Real xi = p(0);
166 
167  switch (deriv_type)
168  {
169  case 0:
170  switch (basis_num)
171  {
172  case 0:
173  return -6*xi + 6*xi*xi;
174  case 1:
175  return 6*xi - 6*xi*xi;
176  case 2:
177  return 1 - 4*xi + 3*xi*xi;
178  case 3:
179  return -2*xi + 3*xi*xi;
180 
181  default:
182  libmesh_error_msg("Invalid shape function index i = " <<
183  basis_num);
184  }
185 
186  default:
187  libmesh_error_msg("Invalid shape function derivative j = " <<
188  deriv_type);
189  }
190 }
191 
192 Real clough_raw_shape(const unsigned int basis_num,
193  const Point & p)
194 {
195  Real xi = p(0);
196 
197  switch (basis_num)
198  {
199  case 0:
200  return 1 - 3*xi*xi + 2*xi*xi*xi;
201  case 1:
202  return 3*xi*xi - 2*xi*xi*xi;
203  case 2:
204  return xi - 2*xi*xi + xi*xi*xi;
205  case 3:
206  return -xi*xi + xi*xi*xi;
207 
208  default:
209  libmesh_error_msg("Invalid shape function index i = " <<
210  basis_num);
211  }
212 }
213 
214 
215 } // end anonymous namespace
216 
217 
218 namespace libMesh
219 {
220 
221 
223 
224 
225 template <>
226 Real FE<1,CLOUGH>::shape(const Elem * elem,
227  const Order order,
228  const unsigned int i,
229  const Point & p,
230  const bool add_p_level)
231 {
232  libmesh_assert(elem);
233 
234  clough_compute_coefs(elem);
235 
236  const ElemType type = elem->type();
237 
238  const Order totalorder =
239  static_cast<Order>(order + add_p_level * elem->p_level());
240 
241  switch (totalorder)
242  {
243  // 3rd-order C1 cubic element
244  case THIRD:
245  {
246  switch (type)
247  {
248  // C1 functions on the C1 cubic edge
249  case EDGE2:
250  case EDGE3:
251  {
252  libmesh_assert_less (i, 4);
253 
254  switch (i)
255  {
256  case 0:
257  return clough_raw_shape(0, p);
258  case 1:
259  return clough_raw_shape(1, p);
260  case 2:
261  return d1xd1x * clough_raw_shape(2, p);
262  case 3:
263  return d2xd2x * clough_raw_shape(3, p);
264  default:
265  libmesh_error_msg("Invalid shape function index i = " << i);
266  }
267  }
268  default:
269  libmesh_error_msg("ERROR: Unsupported element type = " << Utility::enum_to_string(type));
270  }
271  }
272  // by default throw an error
273  default:
274  libmesh_error_msg("ERROR: Unsupported polynomial order = " << totalorder);
275  }
276 }
277 
278 
279 
280 template <>
282  const Order,
283  const unsigned int,
284  const Point &)
285 {
286  libmesh_error_msg("Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
287  return 0.;
288 }
289 
290 template <>
292  const Elem * elem,
293  const unsigned int i,
294  const Point & p,
295  const bool add_p_level)
296 {
297  return FE<1,CLOUGH>::shape(elem, fet.order, i, p, add_p_level);
298 }
299 
300 
301 
302 
303 template <>
305  const Order,
306  const unsigned int,
307  const unsigned int,
308  const Point &)
309 {
310  libmesh_error_msg("Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
311  return 0.;
312 }
313 
314 
315 
316 template <>
318  const Order order,
319  const unsigned int i,
320  const unsigned int j,
321  const Point & p,
322  const bool add_p_level)
323 {
324  libmesh_assert(elem);
325 
326  clough_compute_coefs(elem);
327 
328  const ElemType type = elem->type();
329 
330  const Order totalorder =
331  static_cast<Order>(order + add_p_level * elem->p_level());
332 
333  switch (totalorder)
334  {
335  // 3rd-order C1 cubic element
336  case THIRD:
337  {
338  switch (type)
339  {
340  // C1 functions on the C1 cubic edge
341  case EDGE2:
342  case EDGE3:
343  {
344  switch (i)
345  {
346  case 0:
347  return clough_raw_shape_deriv(0, j, p);
348  case 1:
349  return clough_raw_shape_deriv(1, j, p);
350  case 2:
351  return d1xd1x * clough_raw_shape_deriv(2, j, p);
352  case 3:
353  return d2xd2x * clough_raw_shape_deriv(3, j, p);
354  default:
355  libmesh_error_msg("Invalid shape function index i = " << i);
356  }
357  }
358  default:
359  libmesh_error_msg("ERROR: Unsupported element type = " << Utility::enum_to_string(type));
360  }
361  }
362  // by default throw an error
363  default:
364  libmesh_error_msg("ERROR: Unsupported polynomial order = " << totalorder);
365  }
366 }
367 
368 
369 template <>
371  const Elem * elem,
372  const unsigned int i,
373  const unsigned int j,
374  const Point & p,
375  const bool add_p_level)
376 {
377  return FE<1,CLOUGH>::shape_deriv(elem, fet.order, i, j, p, add_p_level);
378 }
379 
380 
381 
382 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
383 
384 template <>
386  const Order order,
387  const unsigned int i,
388  const unsigned int j,
389  const Point & p,
390  const bool add_p_level)
391 {
392  libmesh_assert(elem);
393 
394  clough_compute_coefs(elem);
395 
396  const ElemType type = elem->type();
397 
398  const Order totalorder =
399  static_cast<Order>(order + add_p_level * elem->p_level());
400 
401  switch (totalorder)
402  {
403  // 3rd-order C1 cubic element
404  case THIRD:
405  {
406  switch (type)
407  {
408  // C1 functions on the C1 cubic edge
409  case EDGE2:
410  case EDGE3:
411  {
412  switch (i)
413  {
414  case 0:
415  return clough_raw_shape_second_deriv(0, j, p);
416  case 1:
417  return clough_raw_shape_second_deriv(1, j, p);
418  case 2:
419  return d1xd1x * clough_raw_shape_second_deriv(2, j, p);
420  case 3:
421  return d2xd2x * clough_raw_shape_second_deriv(3, j, p);
422  default:
423  libmesh_error_msg("Invalid shape function index i = " << i);
424  }
425  }
426  default:
427  libmesh_error_msg("ERROR: Unsupported element type = " << Utility::enum_to_string(type));
428  }
429  }
430  // by default throw an error
431  default:
432  libmesh_error_msg("ERROR: Unsupported polynomial order = " << totalorder);
433  }
434 }
435 
436 
437 template <>
439  const Order,
440  const unsigned int,
441  const unsigned int,
442  const Point &)
443 {
444  libmesh_error_msg("Clough-Tocher elements require the real element \nto construct gradient-based degrees of freedom.");
445  return 0.;
446 }
447 
448 template <>
450  const Elem * elem,
451  const unsigned int i,
452  const unsigned int j,
453  const Point & p,
454  const bool add_p_level)
455 {
456  return FE<1,CLOUGH>::shape_second_deriv(elem, fet.order, i, j, p, add_p_level);
457 }
458 
459 
460 
461 #endif
462 
463 } // namespace libMesh
Real(* shape_deriv_ptr)(const FEType fet, const Elem *elem, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
Typedef for pointer to a function that returns FE shape function derivative values.
Definition: fe_interface.h:623
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
unsigned int n_threads()
Definition: libmesh_base.h:96
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.
Definition: elem.h:94
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
Definition: elem.h:2945
OrderWrapper order
The approximation order of the element.
Definition: fe_type.h:201
The libMesh namespace provides an interface to certain functionality in the library.
LIBMESH_DEFAULT_VECTORIZED_FE(template<>Real FE< 0, BERNSTEIN)
dof_id_type id() const
Definition: dof_object.h:823
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:515
libmesh_assert(ctx)
static const dof_id_type invalid_id
An invalid id to distinguish an uninitialized DofObject.
Definition: dof_object.h:477
static shape_deriv_ptr shape_deriv_function(const unsigned int dim, const FEType &fe_t, const ElemType t)
std::string enum_to_string(const T e)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
FEFamily
defines an enum for finite element families.
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
uint8_t dof_id_type
Definition: id_types.h:67
static FEFamily map_fe_type(const Elem &elem)
Definition: fe_map.C:45