libMesh
fe_hermite_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 #include "libmesh/utility.h"
26 
27 namespace
28 {
29 using namespace libMesh;
30 
31 // Compute the static coefficients for an element
32 void hermite_compute_coefs(const Elem * elem, Real & d1xd1x, Real & d2xd2x)
33 {
34  const FEFamily mapping_family = FEMap::map_fe_type(*elem);
35  const Order mapping_order (elem->default_order());
36  const ElemType mapping_elem_type (elem->type());
37 
38  const FEType map_fe_type(mapping_order, mapping_family);
39 
40  const int n_mapping_shape_functions =
41  FEInterface::n_shape_functions (1, map_fe_type,
42  mapping_elem_type);
43 
44  // Degrees of freedom are at vertices and edge midpoints
45  std::vector<Point> dofpt;
46  dofpt.push_back(Point(-1));
47  dofpt.push_back(Point(1));
48 
49  // Mapping functions - first derivatives at each dofpt
50  std::vector<Real> dxdxi(2);
51  std::vector<Real> dxidx(2);
52 
53  FEInterface::shape_deriv_ptr shape_deriv_ptr =
54  FEInterface::shape_deriv_function(1, map_fe_type);
55 
56  for (int p = 0; p != 2; ++p)
57  {
58  dxdxi[p] = 0;
59  for (int i = 0; i != n_mapping_shape_functions; ++i)
60  {
61  const Real ddxi = shape_deriv_ptr
62  (elem, mapping_order, i, 0, dofpt[p], false);
63  dxdxi[p] += elem->point(i)(0) * ddxi;
64  }
65  }
66 
67  // Calculate derivative scaling factors
68 
69  d1xd1x = dxdxi[0];
70  d2xd2x = dxdxi[1];
71 }
72 
73 
74 } // end anonymous namespace
75 
76 
77 namespace libMesh
78 {
79 
80 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
81 
82 template<>
83 Real FEHermite<1>::hermite_raw_shape_second_deriv (const unsigned int i, const Real xi)
84 {
85  using Utility::pow;
86 
87  switch (i)
88  {
89  case 0:
90  return 1.5 * xi;
91  case 1:
92  return -1.5 * xi;
93  case 2:
94  return 0.5 * (-1. + 3.*xi);
95  case 3:
96  return 0.5 * (1. + 3.*xi);
97  case 4:
98  return (8.*xi*xi + 4.*(xi*xi-1.))/24.;
99  case 5:
100  return (8.*xi*xi*xi + 12.*xi*(xi*xi-1.))/120.;
101  // case 6:
102  // return (8.*pow<4>(xi) + 20.*xi*xi*(xi*xi-1.) +
103  // 2.*(xi*xi-1)*(xi*xi-1))/720.;
104  default:
105  Real denominator = 720., xipower = 1.;
106  for (unsigned n=6; n != i; ++n)
107  {
108  xipower *= xi;
109  denominator *= (n+1);
110  }
111  return (8.*pow<4>(xi)*xipower +
112  (8.*(i-4)+4.)*xi*xi*xipower*(xi*xi-1.) +
113  (i-4)*(i-5)*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator;
114  }
115 }
116 
117 #endif
118 
119 
120 template<>
121 Real FEHermite<1>::hermite_raw_shape_deriv(const unsigned int i, const Real xi)
122 {
123  switch (i)
124  {
125  case 0:
126  return 0.75 * (-1. + xi*xi);
127  case 1:
128  return 0.75 * (1. - xi*xi);
129  case 2:
130  return 0.25 * (-1. - 2.*xi + 3.*xi*xi);
131  case 3:
132  return 0.25 * (-1. + 2.*xi + 3.*xi*xi);
133  case 4:
134  return 4.*xi * (xi*xi-1.)/24.;
135  case 5:
136  return (4*xi*xi*(xi*xi-1.) + (xi*xi-1.)*(xi*xi-1.))/120.;
137  // case 6:
138  // return (4*xi*xi*xi*(xi*xi-1.) + 2*xi*(xi*xi-1.)*(xi*xi-1.))/720.;
139  default:
140  Real denominator = 720., xipower = 1.;
141  for (unsigned n=6; n != i; ++n)
142  {
143  xipower *= xi;
144  denominator *= (n+1);
145  }
146  return (4*xi*xi*xi*xipower*(xi*xi-1.) +
147  (i-4)*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator;
148  }
149 }
150 
151 template<>
152 Real FEHermite<1>::hermite_raw_shape(const unsigned int i, const Real xi)
153 {
154  switch (i)
155  {
156  case 0:
157  return 0.25 * (2. - 3.*xi + xi*xi*xi);
158  case 1:
159  return 0.25 * (2. + 3.*xi - xi*xi*xi);
160  case 2:
161  return 0.25 * (1. - xi - xi*xi + xi*xi*xi);
162  case 3:
163  return 0.25 * (-1. - xi + xi*xi + xi*xi*xi);
164  // All high order terms have the form x^(p-4)(x^2-1)^2/p!
165  case 4:
166  return (xi*xi-1.) * (xi*xi-1.)/24.;
167  case 5:
168  return xi * (xi*xi-1.) * (xi*xi-1.)/120.;
169  // case 6:
170  // return xi*xi * (xi*xi-1.) * (xi*xi-1.)/720.;
171  default:
172  Real denominator = 720., xipower = 1.;
173  for (unsigned n=6; n != i; ++n)
174  {
175  xipower *= xi;
176  denominator *= (n+1);
177  }
178  return (xi*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator;
179 
180  }
181 }
182 
183 
184 template <>
186  const Order,
187  const unsigned int,
188  const Point &)
189 {
190  libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
191  return 0.;
192 }
193 
194 
195 
196 template <>
198  const Order libmesh_dbg_var(order),
199  const unsigned int i,
200  const Point & p,
201  const bool libmesh_dbg_var(add_p_level))
202 {
203  libmesh_assert(elem);
204 
205  // Coefficient naming: d(1)d(2n) is the coefficient of the
206  // global shape function corresponding to value 1 in terms of the
207  // local shape function corresponding to normal derivative 2
208  Real d1xd1x, d2xd2x;
209 
210  hermite_compute_coefs(elem, d1xd1x, d2xd2x);
211 
212  const ElemType type = elem->type();
213 
214 #ifndef NDEBUG
215  const unsigned int totalorder =
216  order + add_p_level * elem->p_level();
217 #endif
218 
219  switch (type)
220  {
221  // C1 functions on the C1 cubic edge
222  case EDGE2:
223  case EDGE3:
224  {
225  libmesh_assert_less (i, totalorder+1);
226 
227  switch (i)
228  {
229  case 0:
230  return FEHermite<1>::hermite_raw_shape(0, p(0));
231  case 1:
232  return d1xd1x * FEHermite<1>::hermite_raw_shape(2, p(0));
233  case 2:
234  return FEHermite<1>::hermite_raw_shape(1, p(0));
235  case 3:
236  return d2xd2x * FEHermite<1>::hermite_raw_shape(3, p(0));
237  default:
238  return FEHermite<1>::hermite_raw_shape(i, p(0));
239  }
240  }
241  default:
242  libmesh_error_msg("ERROR: Unsupported element type = " << type);
243  }
244 }
245 
246 
247 
248 template <>
250  const Order,
251  const unsigned int,
252  const unsigned int,
253  const Point &)
254 {
255  libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
256  return 0.;
257 }
258 
259 
260 
261 template <>
263  const Order libmesh_dbg_var(order),
264  const unsigned int i,
265  const unsigned int,
266  const Point & p,
267  const bool libmesh_dbg_var(add_p_level))
268 {
269  libmesh_assert(elem);
270 
271  // Coefficient naming: d(1)d(2n) is the coefficient of the
272  // global shape function corresponding to value 1 in terms of the
273  // local shape function corresponding to normal derivative 2
274  Real d1xd1x, d2xd2x;
275 
276  hermite_compute_coefs(elem, d1xd1x, d2xd2x);
277 
278  const ElemType type = elem->type();
279 
280 #ifndef NDEBUG
281  const unsigned int totalorder =
282  order + add_p_level * elem->p_level();
283 #endif
284 
285  switch (type)
286  {
287  // C1 functions on the C1 cubic edge
288  case EDGE2:
289  case EDGE3:
290  {
291  libmesh_assert_less (i, totalorder+1);
292 
293  switch (i)
294  {
295  case 0:
297  case 1:
298  return d1xd1x * FEHermite<1>::hermite_raw_shape_deriv(2, p(0));
299  case 2:
301  case 3:
302  return d2xd2x * FEHermite<1>::hermite_raw_shape_deriv(3, p(0));
303  default:
305  }
306  }
307  default:
308  libmesh_error_msg("ERROR: Unsupported element type = " << type);
309  }
310 }
311 
312 
313 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
314 
315 template <>
317  const Order,
318  const unsigned int,
319  const unsigned int,
320  const Point &)
321 {
322  libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
323  return 0.;
324 }
325 
326 
327 template <>
329  const Order libmesh_dbg_var(order),
330  const unsigned int i,
331  const unsigned int,
332  const Point & p,
333  const bool libmesh_dbg_var(add_p_level))
334 {
335  libmesh_assert(elem);
336 
337  // Coefficient naming: d(1)d(2n) is the coefficient of the
338  // global shape function corresponding to value 1 in terms of the
339  // local shape function corresponding to normal derivative 2
340  Real d1xd1x, d2xd2x;
341 
342  hermite_compute_coefs(elem, d1xd1x, d2xd2x);
343 
344  const ElemType type = elem->type();
345 
346 #ifndef NDEBUG
347  const unsigned int totalorder =
348  order + add_p_level * elem->p_level();
349 #endif
350 
351  switch (type)
352  {
353  // C1 functions on the C1 cubic edge
354  case EDGE2:
355  case EDGE3:
356  {
357  libmesh_assert_less (i, totalorder+1);
358 
359  switch (i)
360  {
361  case 0:
363  case 1:
364  return d1xd1x * FEHermite<1>::hermite_raw_shape_second_deriv(2, p(0));
365  case 2:
367  case 3:
368  return d2xd2x * FEHermite<1>::hermite_raw_shape_second_deriv(3, p(0));
369  default:
371  }
372  }
373  default:
374  libmesh_error_msg("ERROR: Unsupported element type = " << type);
375  }
376 }
377 
378 #endif
379 
380 } // 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::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
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::FEHermite::hermite_raw_shape_second_deriv
static Real hermite_raw_shape_second_deriv(const unsigned int basis_num, const Real xi)
1D hermite functions on unit interval
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
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::FEHermite::hermite_raw_shape
static Real hermite_raw_shape(const unsigned int basis_num, const Real xi)
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::Elem::default_order
virtual Order default_order() const =0
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::Utility::pow
T pow(const T &x)
Definition: utility.h:246
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::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::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::FEHermite::hermite_raw_shape_deriv
static Real hermite_raw_shape_deriv(const unsigned int basis_num, const Real xi)
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33