libMesh
fe_lagrange_shape_1D.h
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 // Local includes
20 #include "libmesh/enum_elem_type.h"
21 #include "libmesh/enum_order.h"
22 #include "libmesh/point.h"
23 
24 // Inline functions useful to inline on tensor elements.
25 
26 namespace libMesh
27 {
28 
29 inline
30 Real fe_lagrange_1D_linear_shape(const unsigned int i,
31  const Real xi)
32 {
33  libmesh_assert_less (i, 2);
34 
35  switch (i)
36  {
37  case 0:
38  return .5*(1. - xi);
39 
40  case 1:
41  return .5*(1. + xi);
42 
43  default:
44  libmesh_error_msg("Invalid shape function index i = " << i);
45  }
46 }
47 
48 
49 
50 inline
51 Real fe_lagrange_1D_quadratic_shape(const unsigned int i,
52  const Real xi)
53 {
54  libmesh_assert_less (i, 3);
55 
56  switch (i)
57  {
58  case 0:
59  return .5*xi*(xi - 1.);
60 
61  case 1:
62  return .5*xi*(xi + 1);
63 
64  case 2:
65  return (1. - xi*xi);
66 
67  default:
68  libmesh_error_msg("Invalid shape function index i = " << i);
69  }
70 }
71 
72 
73 
74 inline
75 Real fe_lagrange_1D_cubic_shape(const unsigned int i,
76  const Real xi)
77 {
78  libmesh_assert_less (i, 4);
79 
80  switch (i)
81  {
82  case 0:
83  return 9./16.*(1./9.-xi*xi)*(xi-1.);
84 
85  case 1:
86  return -9./16.*(1./9.-xi*xi)*(xi+1.);
87 
88  case 2:
89  return 27./16.*(1.-xi*xi)*(1./3.-xi);
90 
91  case 3:
92  return 27./16.*(1.-xi*xi)*(1./3.+xi);
93 
94  default:
95  libmesh_error_msg("Invalid shape function index i = " << i);
96  }
97 }
98 
99 
100 
101 inline
103  const unsigned int i,
104  const Real xi)
105 {
106  switch (order)
107  {
108  // Lagrange linears
109  case FIRST:
110  return fe_lagrange_1D_linear_shape(i, xi);
111 
112  // Lagrange quadratics
113  case SECOND:
114  return fe_lagrange_1D_quadratic_shape(i, xi);
115 
116  // Lagrange cubics
117  case THIRD:
118  return fe_lagrange_1D_cubic_shape(i, xi);
119 
120  default:
121  libmesh_error_msg("ERROR: Unsupported polynomial order = " << order);
122  }
123 }
124 
125 
126 
127 inline
129  const unsigned int libmesh_dbg_var(j),
130  const Real)
131 {
132  // only d()/dxi in 1D!
133  libmesh_assert_equal_to (j, 0);
134 
135  libmesh_assert_less (i, 2);
136 
137  switch (i)
138  {
139  case 0:
140  return -.5;
141 
142  case 1:
143  return .5;
144 
145  default:
146  libmesh_error_msg("Invalid shape function index i = " << i);
147  }
148 }
149 
150 
151 inline
153  const unsigned int libmesh_dbg_var(j),
154  const Real xi)
155 {
156  // only d()/dxi in 1D!
157  libmesh_assert_equal_to (j, 0);
158 
159  libmesh_assert_less (i, 3);
160 
161  switch (i)
162  {
163  case 0:
164  return xi-.5;
165 
166  case 1:
167  return xi+.5;
168 
169  case 2:
170  return -2.*xi;
171 
172  default:
173  libmesh_error_msg("Invalid shape function index i = " << i);
174  }
175 }
176 
177 
178 inline
180  const unsigned int libmesh_dbg_var(j),
181  const Real xi)
182 {
183  // only d()/dxi in 1D!
184  libmesh_assert_equal_to (j, 0);
185 
186  libmesh_assert_less (i, 4);
187 
188  switch (i)
189  {
190  case 0:
191  return -9./16.*(3.*xi*xi-2.*xi-1./9.);
192 
193  case 1:
194  return -9./16.*(-3.*xi*xi-2.*xi+1./9.);
195 
196  case 2:
197  return 27./16.*(3.*xi*xi-2./3.*xi-1.);
198 
199  case 3:
200  return 27./16.*(-3.*xi*xi-2./3.*xi+1.);
201 
202  default:
203  libmesh_error_msg("Invalid shape function index i = " << i);
204  }
205 }
206 
207 
208 
209 inline
211  const unsigned int i,
212  const unsigned int j,
213  const Real xi)
214 {
215  switch (order)
216  {
217  case FIRST:
218  return fe_lagrange_1D_linear_shape_deriv(i, j, xi);
219 
220  case SECOND:
221  return fe_lagrange_1D_quadratic_shape_deriv(i, j, xi);
222 
223  case THIRD:
224  return fe_lagrange_1D_cubic_shape_deriv(i, j, xi);
225 
226  default:
227  libmesh_error_msg("ERROR: Unsupported polynomial order = " << order);
228  }
229 }
230 
231 
232 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
233 
234 // fe_lagrange_1D_linear_shape_second_deriv is 0
235 
236 
237 inline
239  const unsigned int libmesh_dbg_var(j),
240  const Real)
241 {
242  // Don't need to switch on j. 1D shape functions
243  // depend on xi only!
244  libmesh_assert_equal_to (j, 0);
245 
246  switch (i)
247  {
248  case 0:
249  return 1.;
250 
251  case 1:
252  return 1.;
253 
254  case 2:
255  return -2.;
256 
257  default:
258  libmesh_error_msg("Invalid shape function index i = " << i);
259  }
260 }
261 
262 
263 inline
265  const unsigned int libmesh_dbg_var(j),
266  const Real xi)
267 {
268  // Don't need to switch on j. 1D shape functions
269  // depend on xi only!
270  libmesh_assert_equal_to (j, 0);
271 
272  switch (i)
273  {
274  case 0:
275  return -9./16.*(6.*xi-2);
276 
277  case 1:
278  return -9./16.*(-6*xi-2.);
279 
280  case 2:
281  return 27./16.*(6*xi-2./3.);
282 
283  case 3:
284  return 27./16.*(-6*xi-2./3.);
285 
286  default:
287  libmesh_error_msg("Invalid shape function index i = " << i);
288  }
289 }
290 
291 
292 
293 inline
295  const unsigned int i,
296  const unsigned int j,
297  const Real xi)
298 {
299  switch (order)
300  {
301  // All second derivatives of linears are zero....
302  case FIRST:
303  return 0.;
304 
305  case SECOND:
307 
308  case THIRD:
310 
311  default:
312  libmesh_error_msg("ERROR: Unsupported polynomial order = " << order);
313  } // end switch (order)
314 }
315 
316 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
317 
318 }
libMesh::fe_lagrange_1D_linear_shape
Real fe_lagrange_1D_linear_shape(const unsigned int i, const Real xi)
Definition: fe_lagrange_shape_1D.h:30
libMesh::fe_lagrange_1D_shape_second_deriv
Real fe_lagrange_1D_shape_second_deriv(const Order order, const unsigned int i, const unsigned int j, const Real xi)
Definition: fe_lagrange_shape_1D.h:294
libMesh::fe_lagrange_1D_cubic_shape_second_deriv
Real fe_lagrange_1D_cubic_shape_second_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real xi)
Definition: fe_lagrange_shape_1D.h:264
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::fe_lagrange_1D_quadratic_shape
Real fe_lagrange_1D_quadratic_shape(const unsigned int i, const Real xi)
Definition: fe_lagrange_shape_1D.h:51
libMesh::fe_lagrange_1D_shape_deriv
Real fe_lagrange_1D_shape_deriv(const Order order, const unsigned int i, const unsigned int j, const Real xi)
Definition: fe_lagrange_shape_1D.h:210
libMesh::SECOND
Definition: enum_order.h:43
libMesh::fe_lagrange_1D_shape
Real fe_lagrange_1D_shape(const Order order, const unsigned int i, const Real xi)
Definition: fe_lagrange_shape_1D.h:102
libMesh::fe_lagrange_1D_cubic_shape
Real fe_lagrange_1D_cubic_shape(const unsigned int i, const Real xi)
Definition: fe_lagrange_shape_1D.h:75
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::FIRST
Definition: enum_order.h:42
libMesh::fe_lagrange_1D_quadratic_shape_deriv
Real fe_lagrange_1D_quadratic_shape_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real xi)
Definition: fe_lagrange_shape_1D.h:152
libMesh::fe_lagrange_1D_quadratic_shape_second_deriv
Real fe_lagrange_1D_quadratic_shape_second_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real)
Definition: fe_lagrange_shape_1D.h:238
libMesh::fe_lagrange_1D_linear_shape_deriv
Real fe_lagrange_1D_linear_shape_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real)
Definition: fe_lagrange_shape_1D.h:128
libMesh::fe_lagrange_1D_cubic_shape_deriv
Real fe_lagrange_1D_cubic_shape_deriv(const unsigned int i, const unsigned int libmesh_dbg_var(j), const Real xi)
Definition: fe_lagrange_shape_1D.h:179