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