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 : // C++ includes
20 :
21 : // Local includes
22 : #include "libmesh/fe.h"
23 : #include "libmesh/elem.h"
24 :
25 :
26 : namespace libMesh
27 : {
28 :
29 :
30 35916 : LIBMESH_DEFAULT_VECTORIZED_FE(1,MONOMIAL)
31 :
32 :
33 : template <>
34 28380 : Real FE<1,MONOMIAL>::shape(const ElemType,
35 : const Order libmesh_dbg_var(order),
36 : const unsigned int i,
37 : const Point & p)
38 : {
39 28380 : const Real xi = p(0);
40 :
41 2317 : libmesh_assert_less_equal (i, static_cast<unsigned int>(order));
42 :
43 : // monomials. since they are hierarchic we only need one case block.
44 28380 : switch (i)
45 : {
46 654 : case 0:
47 654 : return 1.;
48 :
49 6156 : case 1:
50 6156 : return xi;
51 :
52 5208 : case 2:
53 5208 : return xi*xi;
54 :
55 4140 : case 3:
56 4140 : return xi*xi*xi;
57 :
58 2916 : case 4:
59 2916 : return xi*xi*xi*xi;
60 :
61 128 : default:
62 128 : Real val = 1.;
63 9216 : for (unsigned int index = 0; index != i; ++index)
64 7936 : val *= xi;
65 128 : return val;
66 : }
67 : }
68 :
69 :
70 :
71 : template <>
72 28380 : Real FE<1,MONOMIAL>::shape(const Elem * elem,
73 : const Order order,
74 : const unsigned int i,
75 : const Point & p,
76 : const bool add_p_level)
77 : {
78 2317 : libmesh_assert(elem);
79 :
80 33014 : return FE<1,MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i, p);
81 : }
82 :
83 :
84 : template <>
85 0 : Real FE<1,MONOMIAL>::shape(const FEType fet,
86 : const Elem * elem,
87 : const unsigned int i,
88 : const Point & p,
89 : const bool add_p_level)
90 : {
91 0 : libmesh_assert(elem);
92 0 : return FE<1,MONOMIAL>::shape(elem->type(), fet.order + add_p_level*elem->p_level(), i, p);
93 : }
94 :
95 :
96 : template <>
97 15240 : Real FE<1,MONOMIAL>::shape_deriv(const ElemType,
98 : const Order libmesh_dbg_var(order),
99 : const unsigned int i,
100 : const unsigned int libmesh_dbg_var(j),
101 : const Point & p)
102 : {
103 : // only d()/dxi in 1D!
104 :
105 1270 : libmesh_assert_equal_to (j, 0);
106 :
107 15240 : const Real xi = p(0);
108 :
109 1270 : libmesh_assert_less_equal (i, static_cast<unsigned int>(order));
110 :
111 : // monomials. since they are hierarchic we only need one case block.
112 15240 : switch (i)
113 : {
114 310 : case 0:
115 310 : return 0.;
116 :
117 3720 : case 1:
118 3720 : return 1.;
119 :
120 3048 : case 2:
121 3048 : return 2.*xi;
122 :
123 2340 : case 3:
124 2340 : return 3.*xi*xi;
125 :
126 1596 : case 4:
127 1596 : return 4.*xi*xi*xi;
128 :
129 816 : default:
130 816 : Real val = i;
131 4080 : for (unsigned int index = 1; index != i; ++index)
132 3400 : val *= xi;
133 68 : return val;
134 : }
135 : }
136 :
137 :
138 :
139 : template <>
140 15240 : Real FE<1,MONOMIAL>::shape_deriv(const Elem * elem,
141 : const Order order,
142 : const unsigned int i,
143 : const unsigned int j,
144 : const Point & p,
145 : const bool add_p_level)
146 : {
147 1270 : libmesh_assert(elem);
148 :
149 16510 : return FE<1,MONOMIAL>::shape_deriv(elem->type(),
150 16510 : order + add_p_level*elem->p_level(), i, j, p);
151 : }
152 :
153 :
154 : template <>
155 0 : Real FE<1,MONOMIAL>::shape_deriv(const FEType fet,
156 : const Elem * elem,
157 : const unsigned int i,
158 : const unsigned int j,
159 : const Point & p,
160 : const bool add_p_level)
161 : {
162 0 : libmesh_assert(elem);
163 0 : return FE<1,MONOMIAL>::shape_deriv(elem->type(), fet.order + add_p_level*elem->p_level(), i, j, p);
164 :
165 :
166 : }
167 :
168 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
169 :
170 : template <>
171 15240 : Real FE<1,MONOMIAL>::shape_second_deriv(const ElemType,
172 : const Order libmesh_dbg_var(order),
173 : const unsigned int i,
174 : const unsigned int libmesh_dbg_var(j),
175 : const Point & p)
176 : {
177 : // only d()/dxi in 1D!
178 :
179 1270 : libmesh_assert_equal_to (j, 0);
180 :
181 15240 : const Real xi = p(0);
182 :
183 1270 : libmesh_assert_less_equal (i, static_cast<unsigned int>(order));
184 :
185 15240 : switch (i)
186 : {
187 620 : case 0:
188 : case 1:
189 620 : return 0.;
190 :
191 3048 : case 2:
192 3048 : return 2.;
193 :
194 2340 : case 3:
195 2340 : return 6.*xi;
196 :
197 1596 : case 4:
198 1596 : return 12.*xi*xi;
199 :
200 68 : default:
201 68 : Real val = 2.;
202 3264 : for (unsigned int index = 2; index != i; ++index)
203 2584 : val *= (index+1) * xi;
204 68 : return val;
205 : }
206 : }
207 :
208 :
209 :
210 : template <>
211 15240 : Real FE<1,MONOMIAL>::shape_second_deriv(const Elem * elem,
212 : const Order order,
213 : const unsigned int i,
214 : const unsigned int j,
215 : const Point & p,
216 : const bool add_p_level)
217 : {
218 1270 : libmesh_assert(elem);
219 :
220 16510 : return FE<1,MONOMIAL>::shape_second_deriv(elem->type(),
221 16510 : order + add_p_level*elem->p_level(), i, j, p);
222 : }
223 :
224 : template <>
225 0 : Real FE<1,MONOMIAL>::shape_second_deriv(const FEType fet,
226 : const Elem * elem,
227 : const unsigned int i,
228 : const unsigned int j,
229 : const Point & p,
230 : const bool add_p_level)
231 : {
232 0 : libmesh_assert(elem);
233 0 : return FE<1,MONOMIAL>::shape_second_deriv(elem->type(), fet.order + add_p_level*elem->p_level(), i, j, p);
234 : }
235 :
236 : #endif
237 :
238 : } // namespace libMesh
|