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 :
20 : // Local includes
21 : #include "libmesh/elem.h"
22 : #include "libmesh/enum_to_string.h"
23 : #include "libmesh/fe.h"
24 : #include "libmesh/fe_interface.h"
25 : #include "libmesh/fe_macro.h"
26 : #include "libmesh/int_range.h"
27 : #include "libmesh/libmesh_logging.h"
28 :
29 : namespace libMesh
30 : {
31 :
32 : // ------------------------------------------------------------
33 : // XYZ-specific implementations
34 :
35 : // Anonymous namespace for local helper functions
36 : namespace {
37 :
38 47664 : void xyz_nodal_soln(const Elem * elem,
39 : const Order order,
40 : const std::vector<Number> & elem_soln,
41 : std::vector<Number> & nodal_soln,
42 : const bool add_p_level)
43 : {
44 47664 : const unsigned int n_nodes = elem->n_nodes();
45 :
46 47664 : nodal_soln.resize(n_nodes);
47 :
48 51636 : const Order totalorder = order + add_p_level*elem->p_level();
49 :
50 47664 : switch (totalorder)
51 : {
52 : // Constant shape functions
53 0 : case CONSTANT:
54 : {
55 0 : libmesh_assert_equal_to (elem_soln.size(), 1);
56 :
57 0 : std::fill(nodal_soln.begin(), nodal_soln.end(), elem_soln[0]);
58 :
59 0 : return;
60 : }
61 :
62 :
63 : // For other orders do interpolation at the nodes
64 : // explicitly.
65 3972 : default:
66 : {
67 : // FEType object to be passed to various FEInterface functions below.
68 3972 : FEType fe_type(order, XYZ);
69 :
70 : const unsigned int n_sf =
71 47664 : FEInterface::n_shape_functions(fe_type, elem);
72 3972 : libmesh_assert_equal_to (elem_soln.size(), n_sf);
73 :
74 : // Zero before summation
75 3972 : std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
76 :
77 628560 : for (unsigned int n=0; n<n_nodes; n++)
78 : // u_i = Sum (alpha_i phi_i)
79 4485168 : for (unsigned int i=0; i<n_sf; i++)
80 4229628 : nodal_soln[n] += elem_soln[i] *
81 4229628 : FEInterface::shape(fe_type, elem, i, elem->point(n));
82 :
83 3972 : return;
84 : } // default
85 : } // switch
86 : } // xyz_nodal_soln()
87 :
88 :
89 : } // anonymous namespace
90 :
91 :
92 :
93 :
94 :
95 :
96 :
97 : template <unsigned int Dim>
98 1362041 : void FEXYZ<Dim>::init_shape_functions(const std::vector<Point> & qp,
99 : const Elem * elem)
100 : {
101 113269 : libmesh_assert(elem);
102 :
103 : // FIXME: Is this redundant here? Who's calling init_shape_functions
104 : // from code that hasn't already done a determine_calculations()?
105 1362041 : this->determine_calculations();
106 :
107 : // Start logging the shape function initialization
108 226538 : LOG_SCOPE("init_shape_functions()", "FE");
109 :
110 : // The number of quadrature points.
111 226538 : const std::size_t n_qp = qp.size();
112 :
113 : // Number of shape functions in the finite element approximation
114 : // space.
115 226538 : const unsigned int n_approx_shape_functions =
116 1135503 : this->n_dofs(elem,
117 : this->get_order());
118 :
119 : // resize the vectors to hold current data
120 : // Phi are the shape functions used for the FE approximation
121 : {
122 : // (note: GCC 3.4.0 requires the use of this-> here)
123 1362041 : if (this->calculate_phi)
124 1075193 : this->phi.resize (n_approx_shape_functions);
125 1362041 : if (this->calculate_dphi)
126 : {
127 1056401 : this->dphi.resize (n_approx_shape_functions);
128 1056401 : this->dphidx.resize (n_approx_shape_functions);
129 1056401 : this->dphidy.resize (n_approx_shape_functions);
130 1056401 : this->dphidz.resize (n_approx_shape_functions);
131 : }
132 :
133 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
134 1362041 : if (this->calculate_d2phi)
135 : {
136 913817 : this->d2phi.resize (n_approx_shape_functions);
137 913817 : this->d2phidx2.resize (n_approx_shape_functions);
138 913817 : this->d2phidxdy.resize (n_approx_shape_functions);
139 913817 : this->d2phidxdz.resize (n_approx_shape_functions);
140 913817 : this->d2phidy2.resize (n_approx_shape_functions);
141 913817 : this->d2phidydz.resize (n_approx_shape_functions);
142 913817 : this->d2phidz2.resize (n_approx_shape_functions);
143 913817 : this->d2phidxi2.resize (n_approx_shape_functions);
144 : }
145 : #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
146 :
147 19346661 : for (unsigned int i=0; i<n_approx_shape_functions; i++)
148 : {
149 17984620 : if (this->calculate_phi)
150 17457875 : this->phi[i].resize (n_qp);
151 17984620 : if (this->calculate_dphi)
152 : {
153 17120915 : this->dphi[i].resize (n_qp);
154 17120915 : this->dphidx[i].resize (n_qp);
155 17120915 : this->dphidy[i].resize (n_qp);
156 17120915 : this->dphidz[i].resize (n_qp);
157 : }
158 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
159 17984620 : if (this->calculate_d2phi)
160 : {
161 16262975 : this->d2phi[i].resize (n_qp);
162 16262975 : this->d2phidx2[i].resize (n_qp);
163 16262975 : this->d2phidxdy[i].resize (n_qp);
164 16262975 : this->d2phidxdz[i].resize (n_qp);
165 16262975 : this->d2phidy2[i].resize (n_qp);
166 16262975 : this->d2phidydz[i].resize (n_qp);
167 16262975 : this->d2phidz2[i].resize (n_qp);
168 : }
169 : #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
170 : }
171 : }
172 :
173 :
174 :
175 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
176 : //------------------------------------------------------------
177 : // Initialize the data fields, which should only be used for infinite
178 : // elements, to some sensible values, so that using a FE with the
179 : // variational formulation of an InfFE, correct element matrices are
180 : // returned
181 :
182 : {
183 326209 : this->weight.resize (n_qp);
184 326209 : this->dweight.resize (n_qp);
185 326209 : this->dphase.resize (n_qp);
186 :
187 2154487 : for (unsigned int p=0; p<n_qp; p++)
188 : {
189 1828278 : this->weight[p] = 1.;
190 1228810 : this->dweight[p].zero();
191 1228810 : this->dphase[p].zero();
192 : }
193 :
194 : }
195 : #endif // ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
196 :
197 1362041 : if (this->calculate_dual)
198 954 : this->init_dual_shape_functions(n_approx_shape_functions, n_qp);
199 1362041 : }
200 :
201 :
202 :
203 :
204 : template <unsigned int Dim>
205 1362041 : void FEXYZ<Dim>::compute_shape_functions (const Elem * elem,
206 : const std::vector<Point> &)
207 : {
208 113269 : libmesh_assert(elem);
209 :
210 : //-------------------------------------------------------------------------
211 : // Compute the shape function values (and derivatives)
212 : // at the Quadrature points. Note that the actual values
213 : // have already been computed via init_shape_functions
214 :
215 : // Start logging the shape function computation
216 226538 : LOG_SCOPE("compute_shape_functions()", "FE");
217 :
218 326209 : const std::vector<Point> & xyz_qp = this->get_xyz();
219 :
220 : // Compute the value of the derivative shape function i at quadrature point p
221 1362041 : switch (this->dim)
222 : {
223 :
224 3024 : case 1:
225 : {
226 3024 : if (this->calculate_phi)
227 13392 : for (auto i : index_range(this->phi))
228 26592 : for (auto p : index_range(this->phi[i]))
229 17524 : this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]);
230 :
231 3024 : if (this->calculate_dphi)
232 11448 : for (auto i : index_range(this->dphi))
233 19248 : for (auto p : index_range(this->dphi[i]))
234 : {
235 11206 : this->dphi[i][p](0) =
236 11206 : this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
237 :
238 11206 : this->dphi[i][p](1) = this->dphidy[i][p] = 0.;
239 12068 : this->dphi[i][p](2) = this->dphidz[i][p] = 0.;
240 : }
241 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
242 3024 : if (this->calculate_d2phi)
243 11448 : for (auto i : index_range(this->d2phi))
244 19248 : for (auto p : index_range(this->d2phi[i]))
245 : {
246 11206 : this->d2phi[i][p](0,0) =
247 11206 : this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
248 :
249 : #if LIBMESH_DIM>1
250 11206 : this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] =
251 9482 : this->d2phi[i][p](1,0) = 0.;
252 11206 : this->d2phi[i][p](1,1) = this->d2phidy2[i][p] = 0.;
253 : #if LIBMESH_DIM>2
254 11206 : this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] =
255 9482 : this->d2phi[i][p](2,0) = 0.;
256 11206 : this->d2phi[i][p](1,2) = this->d2phidydz[i][p] =
257 9482 : this->d2phi[i][p](2,1) = 0.;
258 12068 : this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = 0.;
259 : #endif
260 : #endif
261 : }
262 : #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
263 :
264 : // All done
265 252 : break;
266 : }
267 :
268 197936 : case 2:
269 : {
270 197936 : if (this->calculate_phi)
271 1055400 : for (auto i : index_range(this->phi))
272 3507230 : for (auto p : index_range(this->phi[i]))
273 2804206 : this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]);
274 :
275 197936 : if (this->calculate_dphi)
276 1030884 : for (auto i : index_range(this->dphi))
277 3143486 : for (auto p : index_range(this->dphi[i]))
278 : {
279 2434018 : this->dphi[i][p](0) =
280 2434018 : this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
281 :
282 2434018 : this->dphi[i][p](1) =
283 2434018 : this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
284 :
285 : #if LIBMESH_DIM == 3
286 2244202 : this->dphi[i][p](2) = // can only assign to the Z component if LIBMESH_DIM==3
287 : #endif
288 2434018 : this->dphidz[i][p] = 0.;
289 : }
290 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
291 197936 : if (this->calculate_d2phi)
292 870516 : for (auto i : index_range(this->d2phi))
293 2501606 : for (auto p : index_range(this->d2phi[i]))
294 : {
295 1877650 : this->d2phi[i][p](0,0) =
296 1877650 : this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
297 :
298 1877650 : this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] =
299 1877650 : this->d2phi[i][p](1,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
300 1877650 : this->d2phi[i][p](1,1) =
301 1877650 : this->d2phidy2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]);
302 : #if LIBMESH_DIM>2
303 1877650 : this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] =
304 1583746 : this->d2phi[i][p](2,0) = 0.;
305 1877650 : this->d2phi[i][p](1,2) = this->d2phidydz[i][p] =
306 1583746 : this->d2phi[i][p](2,1) = 0.;
307 2024602 : this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = 0.;
308 : #endif
309 : }
310 : #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
311 :
312 : // All done
313 16744 : break;
314 : }
315 :
316 1161081 : case 3:
317 : {
318 1161081 : if (this->calculate_phi)
319 16122081 : for (auto i : index_range(this->phi))
320 61003344 : for (auto p : index_range(this->phi[i]))
321 49593855 : this->phi[i][p] = FE<Dim,XYZ>::shape (elem, this->fe_type.order, i, xyz_qp[p]);
322 :
323 1161081 : if (this->calculate_dphi)
324 15818709 : for (auto i : index_range(this->dphi))
325 45041160 : for (auto p : index_range(this->dphi[i]))
326 : {
327 32612943 : this->dphi[i][p](0) =
328 32612943 : this->dphidx[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
329 :
330 32612943 : this->dphi[i][p](1) =
331 32612943 : this->dphidy[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
332 :
333 35081178 : this->dphi[i][p](2) =
334 32612943 : this->dphidz[i][p] = FE<Dim,XYZ>::shape_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]);
335 : }
336 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
337 1161081 : if (this->calculate_d2phi)
338 15044565 : for (auto i : index_range(this->d2phi))
339 41059848 : for (auto p : index_range(this->d2phi[i]))
340 : {
341 29018703 : this->d2phi[i][p](0,0) =
342 29018703 : this->d2phidx2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 0, xyz_qp[p]);
343 :
344 29018703 : this->d2phi[i][p](0,1) = this->d2phidxdy[i][p] =
345 29018703 : this->d2phi[i][p](1,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 1, xyz_qp[p]);
346 29018703 : this->d2phi[i][p](1,1) =
347 29018703 : this->d2phidy2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 2, xyz_qp[p]);
348 29018703 : this->d2phi[i][p](0,2) = this->d2phidxdz[i][p] =
349 29018703 : this->d2phi[i][p](2,0) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 3, xyz_qp[p]);
350 29018703 : this->d2phi[i][p](1,2) = this->d2phidydz[i][p] =
351 29018703 : this->d2phi[i][p](2,1) = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 4, xyz_qp[p]);
352 29018703 : this->d2phi[i][p](2,2) = this->d2phidz2[i][p] = FE<Dim,XYZ>::shape_second_deriv (elem, this->fe_type.order, i, 5, xyz_qp[p]);
353 : }
354 : #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
355 :
356 : // All done
357 96273 : break;
358 : }
359 :
360 0 : default:
361 0 : libmesh_error_msg("ERROR: Invalid dimension " << this->dim);
362 : }
363 1362041 : }
364 :
365 :
366 : // Instantiate (side_) nodal_soln() function for every dimension
367 47664 : LIBMESH_FE_NODAL_SOLN(XYZ, xyz_nodal_soln)
368 0 : LIBMESH_FE_SIDE_NODAL_SOLN(XYZ)
369 :
370 :
371 : // Full specialization of n_dofs() function for every dimension
372 0 : template <> unsigned int FE<0,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
373 0 : template <> unsigned int FE<1,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
374 0 : template <> unsigned int FE<2,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
375 0 : template <> unsigned int FE<3,XYZ>::n_dofs(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
376 :
377 0 : template <> unsigned int FE<0,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
378 3312 : template <> unsigned int FE<1,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
379 464880 : template <> unsigned int FE<2,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
380 2747281 : template <> unsigned int FE<3,XYZ>::n_dofs(const Elem * e, const Order o) { return monomial_n_dofs(e, o); }
381 :
382 : // Full specialization of n_dofs_at_node() function for every dimension.
383 : // XYZ FEMs have no dofs at nodes, only element dofs.
384 0 : template <> unsigned int FE<0,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
385 4554 : template <> unsigned int FE<1,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
386 1302990 : template <> unsigned int FE<2,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
387 10594844 : template <> unsigned int FE<3,XYZ>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
388 :
389 0 : template <> unsigned int FE<0,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
390 3321 : template <> unsigned int FE<1,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
391 207783 : template <> unsigned int FE<2,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
392 2763135 : template <> unsigned int FE<3,XYZ>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
393 :
394 : // Full specialization of n_dofs_per_elem() function for every dimension.
395 0 : template <> unsigned int FE<0,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
396 0 : template <> unsigned int FE<1,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
397 0 : template <> unsigned int FE<2,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
398 0 : template <> unsigned int FE<3,XYZ>::n_dofs_per_elem(const ElemType t, const Order o) { return monomial_n_dofs(t, o); }
399 :
400 0 : template <> unsigned int FE<0,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
401 2628 : template <> unsigned int FE<1,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
402 241155 : template <> unsigned int FE<2,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
403 947058 : template <> unsigned int FE<3,XYZ>::n_dofs_per_elem(const Elem & e, const Order o) { return monomial_n_dofs(&e, o); }
404 :
405 : // Full specialization of get_continuity() function for every dimension.
406 0 : template <> FEContinuity FE<0,XYZ>::get_continuity() const { return DISCONTINUOUS; }
407 13260 : template <> FEContinuity FE<1,XYZ>::get_continuity() const { return DISCONTINUOUS; }
408 53866 : template <> FEContinuity FE<2,XYZ>::get_continuity() const { return DISCONTINUOUS; }
409 90545 : template <> FEContinuity FE<3,XYZ>::get_continuity() const { return DISCONTINUOUS; }
410 :
411 : // Full specialization of is_hierarchic() function for every dimension.
412 : // The XYZ shape functions are hierarchic!
413 0 : template <> bool FE<0,XYZ>::is_hierarchic() const { return true; }
414 48 : template <> bool FE<1,XYZ>::is_hierarchic() const { return true; }
415 262 : template <> bool FE<2,XYZ>::is_hierarchic() const { return true; }
416 644 : template <> bool FE<3,XYZ>::is_hierarchic() const { return true; }
417 :
418 : #ifdef LIBMESH_ENABLE_AMR
419 :
420 : // Full specialization of compute_constraints() function for 2D and
421 : // 3D only. There are no constraints for discontinuous elements, so
422 : // we do nothing.
423 0 : template <> void FE<2,XYZ>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem *) {}
424 0 : template <> void FE<3,XYZ>::compute_constraints (DofConstraints &, DofMap &, const unsigned int, const Elem *) {}
425 :
426 : #endif // #ifdef LIBMESH_ENABLE_AMR
427 :
428 : // Full specialization of shapes_need_reinit() function for every dimension.
429 0 : template <> bool FE<0,XYZ>::shapes_need_reinit() const { return true; }
430 576 : template <> bool FE<1,XYZ>::shapes_need_reinit() const { return true; }
431 228776 : template <> bool FE<2,XYZ>::shapes_need_reinit() const { return true; }
432 800342 : template <> bool FE<3,XYZ>::shapes_need_reinit() const { return true; }
433 :
434 :
435 : // Explicit instantiations for non-static FEXYZ member functions.
436 : // These non-static member functions map more naturally to explicit
437 : // instantiations than the functions above:
438 : //
439 : // 1.) Since they are member functions, they rely on
440 : // private/protected member data, and therefore do not work well
441 : // with the "anonymous function call" model we've used above for
442 : // the specializations.
443 : //
444 : // 2.) There is (IMHO) less chance of the linker calling the
445 : // wrong version of one of these member functions, since there is
446 : // only one FEXYZ.
447 : template LIBMESH_EXPORT void FEXYZ<0>::init_shape_functions(const std::vector<Point> &, const Elem *);
448 : template LIBMESH_EXPORT void FEXYZ<1>::init_shape_functions(const std::vector<Point> &, const Elem *);
449 : template LIBMESH_EXPORT void FEXYZ<2>::init_shape_functions(const std::vector<Point> &, const Elem *);
450 : template LIBMESH_EXPORT void FEXYZ<3>::init_shape_functions(const std::vector<Point> &, const Elem *);
451 :
452 : template LIBMESH_EXPORT void FEXYZ<0>::compute_shape_functions(const Elem *,const std::vector<Point> &);
453 : template LIBMESH_EXPORT void FEXYZ<1>::compute_shape_functions(const Elem *,const std::vector<Point> &);
454 : template LIBMESH_EXPORT void FEXYZ<2>::compute_shape_functions(const Elem *,const std::vector<Point> &);
455 : template LIBMESH_EXPORT void FEXYZ<3>::compute_shape_functions(const Elem *,const std::vector<Point> &);
456 :
457 : } // namespace libMesh
|