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  default:
44  return .5*(1. + xi);
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  default:
66  return (1. - xi*xi);
67  }
68 }
69 
70 
71 
72 inline
73 Real fe_lagrange_1D_cubic_shape(const unsigned int i,
74  const Real xi)
75 {
76  libmesh_assert_less (i, 4);
77 
78  switch (i)
79  {
80  case 0:
81  return 9./16.*(1./9.-xi*xi)*(xi-1.);
82 
83  case 1:
84  return -9./16.*(1./9.-xi*xi)*(xi+1.);
85 
86  case 2:
87  return 27./16.*(1.-xi*xi)*(1./3.-xi);
88 
89  // case 3
90  default:
91  return 27./16.*(1.-xi*xi)*(1./3.+xi);
92  }
93 }
94 
95 
96 
97 inline
99  const unsigned int i,
100  const Real xi)
101 {
102  libmesh_assert_less_equal(order, THIRD);
103 
104  switch (order)
105  {
106  // Lagrange linears
107  case FIRST:
108  return fe_lagrange_1D_linear_shape(i, xi);
109 
110  // Lagrange quadratics
111  case SECOND:
112  return fe_lagrange_1D_quadratic_shape(i, xi);
113 
114  // Lagrange cubics
115  // case THIRD
116  default:
117  return fe_lagrange_1D_cubic_shape(i, xi);
118  }
119 }
120 
121 
122 
123 inline
125  const unsigned int libmesh_dbg_var(j),
126  const Real)
127 {
128  // only d()/dxi in 1D!
129  libmesh_assert_equal_to (j, 0);
130 
131  libmesh_assert_less (i, 2);
132 
133  switch (i)
134  {
135  case 0:
136  return -.5;
137 
138  // case 1
139  default:
140  return .5;
141  }
142 }
143 
144 
145 inline
147  const unsigned int libmesh_dbg_var(j),
148  const Real xi)
149 {
150  // only d()/dxi in 1D!
151  libmesh_assert_equal_to (j, 0);
152 
153  libmesh_assert_less (i, 3);
154 
155  switch (i)
156  {
157  case 0:
158  return xi-.5;
159 
160  case 1:
161  return xi+.5;
162 
163  // case 2
164  default:
165  return -2.*xi;
166  }
167 }
168 
169 
170 inline
172  const unsigned int libmesh_dbg_var(j),
173  const Real xi)
174 {
175  // only d()/dxi in 1D!
176  libmesh_assert_equal_to (j, 0);
177 
178  libmesh_assert_less (i, 4);
179 
180  switch (i)
181  {
182  case 0:
183  return -9./16.*(3.*xi*xi-2.*xi-1./9.);
184 
185  case 1:
186  return -9./16.*(-3.*xi*xi-2.*xi+1./9.);
187 
188  case 2:
189  return 27./16.*(3.*xi*xi-2./3.*xi-1.);
190 
191  // case 3
192  default:
193  return 27./16.*(-3.*xi*xi-2./3.*xi+1.);
194  }
195 }
196 
197 
198 
199 inline
201  const unsigned int i,
202  const unsigned int j,
203  const Real xi)
204 {
205  libmesh_assert_less_equal(order, THIRD);
206 
207  switch (order)
208  {
209  case FIRST:
210  return fe_lagrange_1D_linear_shape_deriv(i, j, xi);
211 
212  case SECOND:
213  return fe_lagrange_1D_quadratic_shape_deriv(i, j, xi);
214 
215  // case THIRD
216  default:
217  return fe_lagrange_1D_cubic_shape_deriv(i, j, xi);
218  }
219 }
220 
221 
222 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
223 
224 // fe_lagrange_1D_linear_shape_second_deriv is 0
225 
226 
227 inline
229  const unsigned int libmesh_dbg_var(j),
230  const Real)
231 {
232  // Don't need to switch on j. 1D shape functions
233  // depend on xi only!
234  libmesh_assert_equal_to (j, 0);
235  libmesh_assert_less(i, 3);
236 
237  switch (i)
238  {
239  case 0:
240  return 1.;
241 
242  case 1:
243  return 1.;
244 
245  // case 2
246  default:
247  return -2.;
248  }
249 }
250 
251 
252 inline
254  const unsigned int libmesh_dbg_var(j),
255  const Real xi)
256 {
257  // Don't need to switch on j. 1D shape functions
258  // depend on xi only!
259  libmesh_assert_equal_to (j, 0);
260  libmesh_assert_less(i, 4);
261 
262  switch (i)
263  {
264  case 0:
265  return -9./16.*(6.*xi-2);
266 
267  case 1:
268  return -9./16.*(-6*xi-2.);
269 
270  case 2:
271  return 27./16.*(6*xi-2./3.);
272 
273  // case 2
274  default:
275  return 27./16.*(-6*xi-2./3.);
276  }
277 }
278 
279 
280 
281 inline
283  const unsigned int i,
284  const unsigned int j,
285  const Real xi)
286 {
287  libmesh_assert_less_equal(order, THIRD);
288 
289  switch (order)
290  {
291  // All second derivatives of linears are zero....
292  case FIRST:
293  return 0.;
294 
295  case SECOND:
297 
298  // case THIRD
299  default:
301  } // end switch (order)
302 }
303 
304 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
305 
306 }
307 
308 #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)