libMesh
fe_hermite_shape_1D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 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  libmesh_error_msg("We'll never get here!");
108  return 0.;
109 }
110 
111 
112 
113 template<>
114 Real FEHermite<1>::hermite_raw_shape_deriv(const unsigned int i, const Real xi)
115 {
116  switch (i)
117  {
118  case 0:
119  return 0.75 * (-1. + xi*xi);
120  case 1:
121  return 0.75 * (1. - xi*xi);
122  case 2:
123  return 0.25 * (-1. - 2.*xi + 3.*xi*xi);
124  case 3:
125  return 0.25 * (-1. + 2.*xi + 3.*xi*xi);
126  case 4:
127  return 4.*xi * (xi*xi-1.)/24.;
128  case 5:
129  return (4*xi*xi*(xi*xi-1.) + (xi*xi-1.)*(xi*xi-1.))/120.;
130  // case 6:
131  // return (4*xi*xi*xi*(xi*xi-1.) + 2*xi*(xi*xi-1.)*(xi*xi-1.))/720.;
132  default:
133  Real denominator = 720., xipower = 1.;
134  for (unsigned n=6; n != i; ++n)
135  {
136  xipower *= xi;
137  denominator *= (n+1);
138  }
139  return (4*xi*xi*xi*xipower*(xi*xi-1.) +
140  (i-4)*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator;
141  }
142 
143  libmesh_error_msg("We'll never get here!");
144  return 0.;
145 }
146 
147 template<>
148 Real FEHermite<1>::hermite_raw_shape(const unsigned int i, const Real xi)
149 {
150  switch (i)
151  {
152  case 0:
153  return 0.25 * (2. - 3.*xi + xi*xi*xi);
154  case 1:
155  return 0.25 * (2. + 3.*xi - xi*xi*xi);
156  case 2:
157  return 0.25 * (1. - xi - xi*xi + xi*xi*xi);
158  case 3:
159  return 0.25 * (-1. - xi + xi*xi + xi*xi*xi);
160  // All high order terms have the form x^(p-4)(x^2-1)^2/p!
161  case 4:
162  return (xi*xi-1.) * (xi*xi-1.)/24.;
163  case 5:
164  return xi * (xi*xi-1.) * (xi*xi-1.)/120.;
165  // case 6:
166  // return xi*xi * (xi*xi-1.) * (xi*xi-1.)/720.;
167  default:
168  Real denominator = 720., xipower = 1.;
169  for (unsigned n=6; n != i; ++n)
170  {
171  xipower *= xi;
172  denominator *= (n+1);
173  }
174  return (xi*xi*xipower*(xi*xi-1.)*(xi*xi-1.))/denominator;
175 
176  }
177 
178  libmesh_error_msg("We'll never get here!");
179  return 0.;
180 }
181 
182 
183 template <>
185  const Order,
186  const unsigned int,
187  const Point &)
188 {
189  libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
190  return 0.;
191 }
192 
193 
194 
195 template <>
197  const Order libmesh_dbg_var(order),
198  const unsigned int i,
199  const Point & p)
200 {
201  libmesh_assert(elem);
202 
203  // Coefficient naming: d(1)d(2n) is the coefficient of the
204  // global shape function corresponding to value 1 in terms of the
205  // local shape function corresponding to normal derivative 2
206  Real d1xd1x, d2xd2x;
207 
208  hermite_compute_coefs(elem, d1xd1x, d2xd2x);
209 
210  const ElemType type = elem->type();
211 
212 #ifndef NDEBUG
213  const unsigned int totalorder = order + elem->p_level();
214 #endif
215 
216  switch (type)
217  {
218  // C1 functions on the C1 cubic edge
219  case EDGE2:
220  case EDGE3:
221  {
222  libmesh_assert_less (i, totalorder+1);
223 
224  switch (i)
225  {
226  case 0:
227  return FEHermite<1>::hermite_raw_shape(0, p(0));
228  case 1:
229  return d1xd1x * FEHermite<1>::hermite_raw_shape(2, p(0));
230  case 2:
231  return FEHermite<1>::hermite_raw_shape(1, p(0));
232  case 3:
233  return d2xd2x * FEHermite<1>::hermite_raw_shape(3, p(0));
234  default:
235  return FEHermite<1>::hermite_raw_shape(i, p(0));
236  }
237  }
238  default:
239  libmesh_error_msg("ERROR: Unsupported element type = " << type);
240  }
241 
242  libmesh_error_msg("We'll never get here!");
243  return 0.;
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 {
268  libmesh_assert(elem);
269 
270  // Coefficient naming: d(1)d(2n) is the coefficient of the
271  // global shape function corresponding to value 1 in terms of the
272  // local shape function corresponding to normal derivative 2
273  Real d1xd1x, d2xd2x;
274 
275  hermite_compute_coefs(elem, d1xd1x, d2xd2x);
276 
277  const ElemType type = elem->type();
278 
279 #ifndef NDEBUG
280  const unsigned int totalorder = order + elem->p_level();
281 #endif
282 
283  switch (type)
284  {
285  // C1 functions on the C1 cubic edge
286  case EDGE2:
287  case EDGE3:
288  {
289  libmesh_assert_less (i, totalorder+1);
290 
291  switch (i)
292  {
293  case 0:
295  case 1:
296  return d1xd1x * FEHermite<1>::hermite_raw_shape_deriv(2, p(0));
297  case 2:
299  case 3:
300  return d2xd2x * FEHermite<1>::hermite_raw_shape_deriv(3, p(0));
301  default:
303  }
304  }
305  default:
306  libmesh_error_msg("ERROR: Unsupported element type = " << type);
307  }
308 
309  libmesh_error_msg("We'll never get here!");
310  return 0.;
311 }
312 
313 
314 
315 template <>
317  const Order libmesh_dbg_var(order),
318  const unsigned int i,
319  const unsigned int,
320  const Point & p)
321 {
322  libmesh_assert(elem);
323 
324  // Coefficient naming: d(1)d(2n) is the coefficient of the
325  // global shape function corresponding to value 1 in terms of the
326  // local shape function corresponding to normal derivative 2
327  Real d1xd1x, d2xd2x;
328 
329  hermite_compute_coefs(elem, d1xd1x, d2xd2x);
330 
331  const ElemType type = elem->type();
332 
333 #ifndef NDEBUG
334  const unsigned int totalorder = order + elem->p_level();
335 #endif
336 
337  switch (type)
338  {
339  // C1 functions on the C1 cubic edge
340  case EDGE2:
341  case EDGE3:
342  {
343  libmesh_assert_less (i, totalorder+1);
344 
345  switch (i)
346  {
347  case 0:
349  case 1:
350  return d1xd1x * FEHermite<1>::hermite_raw_shape_second_deriv(2, p(0));
351  case 2:
353  case 3:
354  return d2xd2x * FEHermite<1>::hermite_raw_shape_second_deriv(3, p(0));
355  default:
357  }
358  }
359  default:
360  libmesh_error_msg("ERROR: Unsupported element type = " << type);
361  }
362 
363  libmesh_error_msg("We'll never get here!");
364  return 0.;
365 }
366 
367 } // namespace libMesh
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
virtual ElemType type() const =0
unsigned int p_level() const
Definition: elem.h:2422
ElemType
Defines an enum for geometric element types.
This is the base class from which all geometric element types are derived.
Definition: elem.h:89
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
The libMesh namespace provides an interface to certain functionality in the library.
static Real hermite_raw_shape_second_deriv(const unsigned int basis_num, const Real xi)
1D hermite functions on unit interval
static Real hermite_raw_shape(const unsigned int basis_num, const Real xi)
libmesh_assert(j)
T pow(const T &x)
Definition: utility.h:195
virtual Order default_order() const =0
PetscErrorCode Vec Mat libmesh_dbg_var(j)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Point & point(const unsigned int i) const
Definition: elem.h:1809
virtual unsigned int n_shape_functions() const libmesh_override
Definition: fe.C:36
Order
defines an enum for polynomial orders.
Definition: enum_order.h:32
static Real hermite_raw_shape_deriv(const unsigned int basis_num, const Real xi)
A Point defines a location in LIBMESH_DIM dimensional Real 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)