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