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