Line data Source code
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 4138319406 : Real fe_lagrange_1D_linear_shape(const unsigned int i,
33 : const Real xi)
34 : {
35 411388430 : libmesh_assert_less (i, 2);
36 :
37 4138319406 : switch (i)
38 : {
39 2069159703 : case 0:
40 2094498634 : return .5*(1. - xi);
41 :
42 2069159703 : case 1:
43 2094498634 : return .5*(1. + xi);
44 :
45 0 : default:
46 0 : libmesh_error_msg("Invalid shape function index i = " << i);
47 : }
48 : }
49 :
50 :
51 :
52 : inline
53 22768542129 : Real fe_lagrange_1D_quadratic_shape(const unsigned int i,
54 : const Real xi)
55 : {
56 1801643358 : libmesh_assert_less (i, 3);
57 :
58 22768542129 : switch (i)
59 : {
60 7589283135 : case 0:
61 7605166659 : return .5*xi*(xi - 1.);
62 :
63 7589283135 : case 1:
64 7589283135 : return .5*xi*(xi + 1);
65 :
66 7589975859 : case 2:
67 7668453368 : return (1. - xi*xi);
68 :
69 0 : default:
70 0 : libmesh_error_msg("Invalid shape function index i = " << i);
71 : }
72 : }
73 :
74 :
75 :
76 : inline
77 24256 : Real fe_lagrange_1D_cubic_shape(const unsigned int i,
78 : const Real xi)
79 : {
80 768 : libmesh_assert_less (i, 4);
81 :
82 24256 : switch (i)
83 : {
84 6064 : case 0:
85 6064 : return 9./16.*(1./9.-xi*xi)*(xi-1.);
86 :
87 6064 : case 1:
88 6064 : return -9./16.*(1./9.-xi*xi)*(xi+1.);
89 :
90 6064 : case 2:
91 6064 : return 27./16.*(1.-xi*xi)*(1./3.-xi);
92 :
93 6064 : case 3:
94 6064 : return 27./16.*(1.-xi*xi)*(1./3.+xi);
95 :
96 0 : default:
97 0 : libmesh_error_msg("Invalid shape function index i = " << i);
98 : }
99 : }
100 :
101 :
102 :
103 : inline
104 78270788 : Real fe_lagrange_1D_shape(const Order order,
105 : const unsigned int i,
106 : const Real xi)
107 : {
108 78270788 : switch (order)
109 : {
110 : // Lagrange linears
111 7544674 : case FIRST:
112 7544674 : return fe_lagrange_1D_linear_shape(i, xi);
113 :
114 : // Lagrange quadratics
115 70701858 : case SECOND:
116 70701858 : return fe_lagrange_1D_quadratic_shape(i, xi);
117 :
118 : // Lagrange cubics
119 24256 : case THIRD:
120 24256 : return fe_lagrange_1D_cubic_shape(i, xi);
121 :
122 0 : default:
123 0 : libmesh_error_msg("ERROR: Unsupported polynomial order = " << order);
124 : }
125 : }
126 :
127 :
128 :
129 : inline
130 1261389220 : Real fe_lagrange_1D_linear_shape_deriv(const unsigned int i,
131 : const unsigned int libmesh_dbg_var(j),
132 : const Real)
133 : {
134 : // only d()/dxi in 1D!
135 119711422 : libmesh_assert_equal_to (j, 0);
136 :
137 119711422 : libmesh_assert_less (i, 2);
138 :
139 1261389220 : switch (i)
140 : {
141 59855711 : case 0:
142 59855711 : return -.5;
143 :
144 630694610 : case 1:
145 630694610 : return .5;
146 :
147 0 : default:
148 0 : libmesh_error_msg("Invalid shape function index i = " << i);
149 : }
150 : }
151 :
152 :
153 : inline
154 10712654301 : Real fe_lagrange_1D_quadratic_shape_deriv(const unsigned int i,
155 : const unsigned int libmesh_dbg_var(j),
156 : const Real xi)
157 : {
158 : // only d()/dxi in 1D!
159 829051189 : libmesh_assert_equal_to (j, 0);
160 :
161 829051189 : libmesh_assert_less (i, 3);
162 :
163 10712654301 : switch (i)
164 : {
165 3570811599 : case 0:
166 3580404670 : return xi-.5;
167 :
168 3570811599 : case 1:
169 3570811599 : return xi+.5;
170 :
171 3571031103 : case 2:
172 3608928430 : return -2.*xi;
173 :
174 0 : default:
175 0 : libmesh_error_msg("Invalid shape function index i = " << i);
176 : }
177 : }
178 :
179 :
180 : inline
181 15920 : Real fe_lagrange_1D_cubic_shape_deriv(const unsigned int i,
182 : const unsigned int libmesh_dbg_var(j),
183 : const Real xi)
184 : {
185 : // only d()/dxi in 1D!
186 512 : libmesh_assert_equal_to (j, 0);
187 :
188 512 : libmesh_assert_less (i, 4);
189 :
190 15920 : switch (i)
191 : {
192 3980 : case 0:
193 3980 : return -9./16.*(3.*xi*xi-2.*xi-1./9.);
194 :
195 3980 : case 1:
196 3980 : return -9./16.*(-3.*xi*xi-2.*xi+1./9.);
197 :
198 3980 : case 2:
199 3980 : return 27./16.*(3.*xi*xi-2./3.*xi-1.);
200 :
201 3980 : case 3:
202 3980 : return 27./16.*(-3.*xi*xi-2./3.*xi+1.);
203 :
204 0 : default:
205 0 : libmesh_error_msg("Invalid shape function index i = " << i);
206 : }
207 : }
208 :
209 :
210 :
211 : inline
212 347507413 : Real fe_lagrange_1D_shape_deriv(const Order order,
213 : const unsigned int i,
214 : const unsigned int j,
215 : const Real xi)
216 : {
217 347507413 : switch (order)
218 : {
219 279486968 : case FIRST:
220 279486968 : return fe_lagrange_1D_linear_shape_deriv(i, j, xi);
221 :
222 68004525 : case SECOND:
223 68004525 : return fe_lagrange_1D_quadratic_shape_deriv(i, j, xi);
224 :
225 15920 : case THIRD:
226 15920 : return fe_lagrange_1D_cubic_shape_deriv(i, j, xi);
227 :
228 0 : default:
229 0 : 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
240 163664424 : Real fe_lagrange_1D_quadratic_shape_second_deriv(const unsigned int i,
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 13920806 : libmesh_assert_equal_to (j, 0);
247 :
248 163664424 : switch (i)
249 : {
250 4535512 : case 0:
251 4535512 : return 1.;
252 :
253 4535512 : case 1:
254 4535512 : return 1.;
255 :
256 54577324 : case 2:
257 54577324 : return -2.;
258 :
259 0 : default:
260 0 : libmesh_error_msg("Invalid shape function index i = " << i);
261 : }
262 : }
263 :
264 :
265 : inline
266 0 : Real fe_lagrange_1D_cubic_shape_second_deriv(const unsigned int i,
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 0 : libmesh_assert_equal_to (j, 0);
273 :
274 0 : switch (i)
275 : {
276 0 : case 0:
277 0 : return -9./16.*(6.*xi-2);
278 :
279 0 : case 1:
280 0 : return -9./16.*(-6*xi-2.);
281 :
282 0 : case 2:
283 0 : return 27./16.*(6*xi-2./3.);
284 :
285 0 : case 3:
286 0 : return 27./16.*(-6*xi-2./3.);
287 :
288 0 : default:
289 0 : libmesh_error_msg("Invalid shape function index i = " << i);
290 : }
291 : }
292 :
293 :
294 :
295 : inline
296 1578960 : Real fe_lagrange_1D_shape_second_deriv(const Order order,
297 : const unsigned int i,
298 : const unsigned int j,
299 : const Real xi)
300 : {
301 1578960 : switch (order)
302 : {
303 : // All second derivatives of linears are zero....
304 12388 : case FIRST:
305 12388 : return 0.;
306 :
307 1009146 : case SECOND:
308 1009146 : return fe_lagrange_1D_quadratic_shape_second_deriv(i, j, xi);
309 :
310 0 : case THIRD:
311 0 : return fe_lagrange_1D_cubic_shape_second_deriv(i, j, xi);
312 :
313 0 : default:
314 0 : 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
|