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