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/dof_map.h"
22 : #include "libmesh/elem.h"
23 : #include "libmesh/enum_to_string.h"
24 : #include "libmesh/fe.h"
25 : #include "libmesh/fe_interface.h"
26 : #include "libmesh/fe_macro.h"
27 : #include "libmesh/tensor_value.h"
28 :
29 :
30 : namespace libMesh
31 : {
32 :
33 :
34 0 : LIBMESH_DEFAULT_VECTORIZED_FE(0,HIERARCHIC_VEC)
35 0 : LIBMESH_DEFAULT_VECTORIZED_FE(1,HIERARCHIC_VEC)
36 139988 : LIBMESH_DEFAULT_VECTORIZED_FE(2,HIERARCHIC_VEC)
37 0 : LIBMESH_DEFAULT_VECTORIZED_FE(3,HIERARCHIC_VEC)
38 0 : LIBMESH_DEFAULT_VECTORIZED_FE(0,L2_HIERARCHIC_VEC)
39 0 : LIBMESH_DEFAULT_VECTORIZED_FE(1,L2_HIERARCHIC_VEC)
40 16940520 : LIBMESH_DEFAULT_VECTORIZED_FE(2,L2_HIERARCHIC_VEC)
41 23081472 : LIBMESH_DEFAULT_VECTORIZED_FE(3,L2_HIERARCHIC_VEC)
42 :
43 :
44 : // ------------------------------------------------------------
45 : // Hierarchic-specific implementations
46 :
47 :
48 : // Anonymous namespace for local helper functions
49 : namespace {
50 55752 : void hierarchic_vec_nodal_soln(const Elem * elem,
51 : const Order order,
52 : const std::vector<Number> & elem_soln,
53 : const int dim,
54 : std::vector<Number> & nodal_soln,
55 : const bool add_p_level)
56 : {
57 55752 : const unsigned int n_nodes = elem->n_nodes();
58 :
59 : // Constant shape functions can't be supported, even for
60 : // L2_HIERARCHIC*, without breaking the "HIERARCHIC is
61 : // hierarchic" guarantee
62 4646 : libmesh_assert(order != CONSTANT);
63 :
64 55752 : nodal_soln.resize(dim*n_nodes);
65 :
66 : // Do interpolation at the nodes explicitly.
67 4646 : FEType fe_type(order, HIERARCHIC);
68 :
69 : const unsigned int n_sf =
70 55752 : FEInterface::n_shape_functions(fe_type, elem, add_p_level);
71 :
72 9292 : std::vector<Point> refspace_nodes;
73 55752 : FEBase::get_refspace_nodes(elem->type(),refspace_nodes);
74 4646 : libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
75 4646 : libmesh_assert_equal_to (elem_soln.size(), n_sf*dim);
76 :
77 : // Zero before summation
78 4646 : std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
79 :
80 590280 : for (unsigned int n=0; n<n_nodes; n++)
81 1861632 : for (int d=0; d != dim; ++d)
82 : // u_i = Sum (alpha_i phi_i); we're here only looking
83 : // at vector components in direction d
84 10696752 : for (unsigned int i=0; i<n_sf; i++)
85 10150452 : nodal_soln[n*dim+d] += elem_soln[i*dim+d] *
86 10150452 : FEInterface::shape(fe_type, elem, i, refspace_nodes[n]);
87 :
88 55752 : }// void hierarchic_vec_nodal_soln
89 :
90 : } // anonymous namespace
91 :
92 :
93 : // Do full-specialization for every dimension, instead
94 : // of explicit instantiation at the end of this file.
95 : // This could be macro-ified so that it fits on one line...
96 0 : LIBMESH_FE_NODAL_SOLN_DIM(HIERARCHIC_VEC, (FE<0, HIERARCHIC>::nodal_soln), 0)
97 0 : LIBMESH_FE_NODAL_SOLN_DIM(HIERARCHIC_VEC, (FE<1, HIERARCHIC>::nodal_soln), 1)
98 :
99 : template <>
100 37320 : void FE<2,HIERARCHIC_VEC>::nodal_soln(const Elem * elem,
101 : const Order order,
102 : const std::vector<Number> & elem_soln,
103 : std::vector<Number> & nodal_soln,
104 : const bool add_p_level,
105 : const unsigned)
106 37320 : { hierarchic_vec_nodal_soln(elem, order, elem_soln, 2 /*dim*/, nodal_soln, add_p_level); }
107 :
108 : template <>
109 18432 : void FE<3,HIERARCHIC_VEC>::nodal_soln(const Elem * elem,
110 : const Order order,
111 : const std::vector<Number> & elem_soln,
112 : std::vector<Number> & nodal_soln,
113 : const bool add_p_level,
114 : const unsigned)
115 18432 : { hierarchic_vec_nodal_soln(elem, order, elem_soln, 3 /*dim*/, nodal_soln, add_p_level); }
116 :
117 0 : LIBMESH_FE_SIDE_NODAL_SOLN(HIERARCHIC_VEC)
118 :
119 0 : LIBMESH_FE_NODAL_SOLN_DIM(L2_HIERARCHIC_VEC, (FE<0, HIERARCHIC_VEC>::nodal_soln), 0)
120 0 : LIBMESH_FE_NODAL_SOLN_DIM(L2_HIERARCHIC_VEC, (FE<1, HIERARCHIC_VEC>::nodal_soln), 1)
121 21600 : LIBMESH_FE_NODAL_SOLN_DIM(L2_HIERARCHIC_VEC, (FE<2, HIERARCHIC_VEC>::nodal_soln), 2)
122 18432 : LIBMESH_FE_NODAL_SOLN_DIM(L2_HIERARCHIC_VEC, (FE<3, HIERARCHIC_VEC>::nodal_soln), 3)
123 :
124 0 : LIBMESH_FE_SIDE_NODAL_SOLN(L2_HIERARCHIC_VEC)
125 :
126 :
127 : // Specialize for shape function routines by leveraging scalar HIERARCHIC elements
128 :
129 : // 0-D
130 0 : template <> RealGradient FE<0,HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
131 : const unsigned int i, const Point & p)
132 : {
133 0 : Real value = FE<0,HIERARCHIC>::shape( type, order, i, p );
134 0 : return libMesh::RealGradient( value );
135 : }
136 0 : template <> RealGradient FE<0,HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
137 : const unsigned int i, const unsigned int j,
138 : const Point & p)
139 : {
140 0 : Real value = FE<0,HIERARCHIC>::shape_deriv( type, order, i, j, p );
141 0 : return libMesh::RealGradient( value );
142 : }
143 :
144 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
145 :
146 0 : template <> RealGradient FE<0,HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
147 : const unsigned int i, const unsigned int j,
148 : const Point & p)
149 : {
150 0 : Real value = FE<0,HIERARCHIC>::shape_second_deriv( type, order, i, j, p );
151 0 : return libMesh::RealGradient( value );
152 : }
153 : #endif
154 :
155 0 : template <> RealGradient FE<0,L2_HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
156 : const unsigned int i, const Point & p)
157 : {
158 0 : return FE<0,HIERARCHIC_VEC>::shape(type, order, i, p);
159 : }
160 0 : template <> RealGradient FE<0,L2_HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
161 : const unsigned int i, const unsigned int j,
162 : const Point & p)
163 : {
164 0 : return FE<0,HIERARCHIC_VEC>::shape_deriv(type, order, i, j, p);
165 : }
166 :
167 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
168 :
169 0 : template <> RealGradient FE<0,L2_HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
170 : const unsigned int i, const unsigned int j,
171 : const Point & p)
172 : {
173 0 : return FE<0,HIERARCHIC_VEC>::shape_second_deriv(type, order, i, j, p);
174 : }
175 : #endif
176 :
177 : // 1-D
178 0 : template <> RealGradient FE<1,HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
179 : const unsigned int i, const Point & p)
180 : {
181 0 : Real value = FE<1,HIERARCHIC>::shape( type, order, i, p );
182 0 : return libMesh::RealGradient( value );
183 : }
184 0 : template <> RealGradient FE<1,HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
185 : const unsigned int i, const unsigned int j,
186 : const Point & p)
187 : {
188 0 : Real value = FE<1,HIERARCHIC>::shape_deriv( type, order, i, j, p );
189 0 : return libMesh::RealGradient( value );
190 : }
191 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
192 0 : template <> RealGradient FE<1,HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
193 : const unsigned int i, const unsigned int j,
194 : const Point & p)
195 : {
196 0 : Real value = FE<1,HIERARCHIC>::shape_second_deriv( type, order, i, j, p );
197 0 : return libMesh::RealGradient( value );
198 : }
199 :
200 : #endif
201 :
202 0 : template <> RealGradient FE<1,L2_HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
203 : const unsigned int i, const Point & p)
204 : {
205 0 : return FE<1,HIERARCHIC_VEC>::shape(type, order, i, p);
206 : }
207 0 : template <> RealGradient FE<1,L2_HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
208 : const unsigned int i, const unsigned int j,
209 : const Point & p)
210 : {
211 0 : return FE<1,HIERARCHIC_VEC>::shape_deriv(type, order, i, j, p);
212 : }
213 :
214 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
215 :
216 0 : template <> RealGradient FE<1,L2_HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
217 : const unsigned int i, const unsigned int j,
218 : const Point & p)
219 : {
220 0 : return FE<1,HIERARCHIC_VEC>::shape_second_deriv(type, order, i, j, p);
221 : }
222 : #endif
223 :
224 : // 2-D
225 0 : template <> RealGradient FE<2,HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
226 : const unsigned int i, const Point & p)
227 : {
228 0 : Real value = FE<2,HIERARCHIC>::shape( type, order, i/2, p );
229 :
230 0 : switch( i%2 )
231 : {
232 0 : case 0:
233 0 : return libMesh::RealGradient( value );
234 :
235 0 : case 1:
236 0 : return libMesh::RealGradient( Real(0), value );
237 :
238 0 : default:
239 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
240 : }
241 :
242 : //dummy
243 : return libMesh::RealGradient();
244 : }
245 0 : template <> RealGradient FE<2,HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
246 : const unsigned int i, const unsigned int j,
247 : const Point & p)
248 : {
249 0 : Real value = FE<2,HIERARCHIC>::shape_deriv( type, order, i/2, j, p );
250 :
251 0 : switch( i%2 )
252 : {
253 0 : case 0:
254 0 : return libMesh::RealGradient( value );
255 :
256 0 : case 1:
257 0 : return libMesh::RealGradient( Real(0), value );
258 :
259 0 : default:
260 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
261 : }
262 :
263 : //dummy
264 : return libMesh::RealGradient();
265 : }
266 :
267 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
268 0 : template <> RealGradient FE<2,HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
269 : const unsigned int i, const unsigned int j,
270 : const Point & p)
271 : {
272 0 : Real value = FE<2,HIERARCHIC>::shape_second_deriv( type, order, i/2, j, p );
273 :
274 0 : switch( i%2 )
275 : {
276 0 : case 0:
277 0 : return libMesh::RealGradient( value );
278 :
279 0 : case 1:
280 0 : return libMesh::RealGradient( Real(0), value );
281 :
282 0 : default:
283 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
284 : }
285 :
286 : //dummy
287 : return libMesh::RealGradient();
288 : }
289 :
290 : #endif
291 :
292 0 : template <> RealGradient FE<2,L2_HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
293 : const unsigned int i, const Point & p)
294 : {
295 0 : return FE<2,HIERARCHIC_VEC>::shape(type, order, i, p);
296 : }
297 :
298 0 : template <> RealGradient FE<2,L2_HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
299 : const unsigned int i, const unsigned int j,
300 : const Point & p)
301 : {
302 0 : return FE<2,HIERARCHIC_VEC>::shape_deriv(type, order, i, j, p);
303 : }
304 :
305 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
306 0 : template <> RealGradient FE<2,L2_HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
307 : const unsigned int i, const unsigned int j,
308 : const Point & p)
309 : {
310 0 : return FE<2,HIERARCHIC_VEC>::shape_second_deriv(type, order, i, j, p);
311 : }
312 :
313 : #endif
314 :
315 : // 3-D
316 0 : template <> RealGradient FE<3,HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
317 : const unsigned int i, const Point & p)
318 : {
319 0 : Real value = FE<3,HIERARCHIC>::shape( type, order, i/3, p );
320 :
321 0 : switch( i%3 )
322 : {
323 0 : case 0:
324 0 : return libMesh::RealGradient( value );
325 :
326 0 : case 1:
327 0 : return libMesh::RealGradient( Real(0), value );
328 :
329 0 : case 2:
330 0 : return libMesh::RealGradient( Real(0), Real(0), value );
331 :
332 0 : default:
333 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
334 : }
335 :
336 : //dummy
337 : return libMesh::RealGradient();
338 : }
339 0 : template <> RealGradient FE<3,HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
340 : const unsigned int i, const unsigned int j,
341 : const Point & p)
342 : {
343 0 : Real value = FE<3,HIERARCHIC>::shape_deriv( type, order, i/3, j, p );
344 :
345 0 : switch( i%3 )
346 : {
347 0 : case 0:
348 0 : return libMesh::RealGradient( value );
349 :
350 0 : case 1:
351 0 : return libMesh::RealGradient( Real(0), value );
352 :
353 0 : case 2:
354 0 : return libMesh::RealGradient( Real(0), Real(0), value );
355 :
356 0 : default:
357 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
358 : }
359 :
360 : //dummy
361 : return libMesh::RealGradient();
362 : }
363 :
364 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
365 :
366 0 : template <> RealGradient FE<3,HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
367 : const unsigned int i, const unsigned int j,
368 : const Point & p)
369 : {
370 0 : Real value = FE<3,HIERARCHIC>::shape_second_deriv( type, order, i/3, j, p );
371 :
372 0 : switch( i%3 )
373 : {
374 0 : case 0:
375 0 : return libMesh::RealGradient( value );
376 :
377 0 : case 1:
378 0 : return libMesh::RealGradient( Real(0), value );
379 :
380 0 : case 2:
381 0 : return libMesh::RealGradient( Real(0), Real(0), value );
382 :
383 0 : default:
384 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
385 : }
386 :
387 : //dummy
388 : return libMesh::RealGradient();
389 : }
390 :
391 : #endif
392 :
393 0 : template <> RealGradient FE<3,L2_HIERARCHIC_VEC>::shape(const ElemType type, const Order order,
394 : const unsigned int i, const Point & p)
395 : {
396 0 : return FE<3,HIERARCHIC_VEC>::shape(type, order, i, p);
397 : }
398 :
399 0 : template <> RealGradient FE<3,L2_HIERARCHIC_VEC>::shape_deriv(const ElemType type, const Order order,
400 : const unsigned int i, const unsigned int j,
401 : const Point & p)
402 : {
403 0 : return FE<3,HIERARCHIC_VEC>::shape_deriv(type, order, i, j, p);
404 : }
405 :
406 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
407 0 : template <> RealGradient FE<3,L2_HIERARCHIC_VEC>::shape_second_deriv(const ElemType type, const Order order,
408 : const unsigned int i, const unsigned int j,
409 : const Point & p)
410 : {
411 0 : return FE<3,HIERARCHIC_VEC>::shape_second_deriv(type, order, i, j, p);
412 : }
413 :
414 : #endif
415 :
416 :
417 : // 0-D
418 0 : template <> RealGradient FE<0,HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
419 : const unsigned int i, const Point & p,
420 : const bool add_p_level)
421 : {
422 0 : const Real value = FE<0,HIERARCHIC>::shape(elem, order, i, p, add_p_level);
423 0 : return libMesh::RealGradient( value );
424 : }
425 0 : template <> RealGradient FE<0,HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
426 : const unsigned int i, const unsigned int j,
427 : const Point & p,
428 : const bool add_p_level)
429 : {
430 0 : const Real value = FE<0,HIERARCHIC>::shape_deriv(elem, order, i, j, p, add_p_level);
431 0 : return libMesh::RealGradient( value );
432 : }
433 :
434 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
435 :
436 0 : template <> RealGradient FE<0,HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
437 : const unsigned int i, const unsigned int j,
438 : const Point & p,
439 : const bool add_p_level)
440 : {
441 0 : Real value = FE<0,HIERARCHIC>::shape_second_deriv(elem, order, i, j, p, add_p_level);
442 0 : return libMesh::RealGradient( value );
443 : }
444 :
445 : #endif
446 :
447 0 : template <> RealGradient FE<0,L2_HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
448 : const unsigned int i, const Point & p,
449 : const bool add_p_level)
450 : {
451 0 : return FE<0,HIERARCHIC_VEC>::shape(elem, order, i, p, add_p_level);
452 : }
453 :
454 0 : template <> RealGradient FE<0,L2_HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
455 : const unsigned int i, const unsigned int j,
456 : const Point & p,
457 : const bool add_p_level)
458 : {
459 0 : return FE<0,HIERARCHIC_VEC>::shape_deriv(elem, order, i, j, p, add_p_level);
460 : }
461 :
462 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
463 :
464 0 : template <> RealGradient FE<0,L2_HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
465 : const unsigned int i, const unsigned int j,
466 : const Point & p,
467 : const bool add_p_level)
468 : {
469 0 : return FE<0,HIERARCHIC_VEC>::shape_second_deriv(elem, order, i, j, p, add_p_level);
470 : }
471 :
472 : #endif
473 :
474 : // 1-D
475 0 : template <> RealGradient FE<1,HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
476 : const unsigned int i, const Point & p,
477 : const bool add_p_level)
478 : {
479 0 : Real value = FE<1,HIERARCHIC>::shape(elem, order, i, p, add_p_level);
480 0 : return libMesh::RealGradient( value );
481 : }
482 0 : template <> RealGradient FE<1,HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
483 : const unsigned int i, const unsigned int j,
484 : const Point & p,
485 : const bool add_p_level)
486 : {
487 0 : Real value = FE<1,HIERARCHIC>::shape_deriv(elem, order, i, j, p, add_p_level);
488 0 : return libMesh::RealGradient( value );
489 : }
490 :
491 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
492 0 : template <> RealGradient FE<1,HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
493 : const unsigned int i, const unsigned int j,
494 : const Point & p,
495 : const bool add_p_level)
496 : {
497 0 : Real value = FE<1,HIERARCHIC>::shape_second_deriv(elem, order, i, j, p, add_p_level);
498 0 : return libMesh::RealGradient( value );
499 : }
500 :
501 : #endif
502 :
503 0 : template <> RealGradient FE<1,L2_HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
504 : const unsigned int i, const Point & p,
505 : const bool add_p_level)
506 : {
507 0 : return FE<1,HIERARCHIC_VEC>::shape(elem, order, i, p, add_p_level);
508 : }
509 :
510 0 : template <> RealGradient FE<1,L2_HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
511 : const unsigned int i, const unsigned int j,
512 : const Point & p,
513 : const bool add_p_level)
514 : {
515 0 : return FE<1,HIERARCHIC_VEC>::shape_deriv(elem, order, i, j, p, add_p_level);
516 : }
517 :
518 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
519 :
520 0 : template <> RealGradient FE<1,L2_HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
521 : const unsigned int i, const unsigned int j,
522 : const Point & p,
523 : const bool add_p_level)
524 : {
525 0 : return FE<1,HIERARCHIC_VEC>::shape_second_deriv(elem, order, i, j, p, add_p_level);
526 : }
527 :
528 : #endif
529 :
530 : // 2-D
531 52697848 : template <> RealGradient FE<2,HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
532 : const unsigned int i, const Point & p,
533 : const bool add_p_level)
534 : {
535 52697848 : const Real value = FE<2,HIERARCHIC>::shape(elem, order, i/2, p, add_p_level);
536 :
537 52697848 : switch( i%2 )
538 : {
539 26348924 : case 0:
540 2179976 : return libMesh::RealGradient( value );
541 :
542 26348924 : case 1:
543 2179976 : return libMesh::RealGradient( Real(0), value );
544 :
545 0 : default:
546 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
547 : }
548 :
549 : //dummy
550 : return libMesh::RealGradient();
551 : }
552 27767032 : template <> RealGradient FE<2,HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
553 : const unsigned int i, const unsigned int j,
554 : const Point & p,
555 : const bool add_p_level)
556 : {
557 27767032 : const Real value = FE<2,HIERARCHIC>::shape_deriv(elem, order, i/2, j, p, add_p_level);
558 :
559 27767032 : switch( i%2 )
560 : {
561 13883516 : case 0:
562 1143112 : return libMesh::RealGradient( value );
563 :
564 13883516 : case 1:
565 1143112 : return libMesh::RealGradient( Real(0), value );
566 :
567 0 : default:
568 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
569 : }
570 :
571 : //dummy
572 : return libMesh::RealGradient();
573 : }
574 :
575 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
576 0 : template <> RealGradient FE<2,HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
577 : const unsigned int i, const unsigned int j,
578 : const Point & p,
579 : const bool add_p_level)
580 : {
581 0 : const Real value = FE<2,HIERARCHIC>::shape_second_deriv(elem, order, i/2, j, p, add_p_level);
582 :
583 0 : switch( i%2 )
584 : {
585 0 : case 0:
586 0 : return libMesh::RealGradient( value );
587 :
588 0 : case 1:
589 0 : return libMesh::RealGradient( Real(0), value );
590 :
591 0 : default:
592 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
593 : }
594 :
595 : //dummy
596 : return libMesh::RealGradient();
597 : }
598 :
599 : #endif
600 :
601 51645120 : template <> RealGradient FE<2,L2_HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
602 : const unsigned int i, const Point & p,
603 : const bool add_p_level)
604 : {
605 51645120 : return FE<2,HIERARCHIC_VEC>::shape(elem, order, i, p, add_p_level);
606 : }
607 27125760 : template <> RealGradient FE<2,L2_HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
608 : const unsigned int i, const unsigned int j,
609 : const Point & p,
610 : const bool add_p_level)
611 : {
612 27125760 : return FE<2,HIERARCHIC_VEC>::shape_deriv(elem, order, i, j, p, add_p_level);
613 : }
614 :
615 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
616 :
617 0 : template <> RealGradient FE<2,L2_HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
618 : const unsigned int i, const unsigned int j,
619 : const Point & p,
620 : const bool add_p_level)
621 : {
622 0 : return FE<2,HIERARCHIC_VEC>::shape_second_deriv(elem, order, i, j, p, add_p_level);
623 : }
624 :
625 : #endif
626 :
627 : // 3-D
628 75423744 : template <> RealGradient FE<3,HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
629 : const unsigned int i, const Point & p,
630 : const bool add_p_level)
631 : {
632 75423744 : const Real value = FE<3,HIERARCHIC>::shape(elem, order, i/3, p, add_p_level);
633 :
634 75423744 : switch( i%3 )
635 : {
636 25141248 : case 0:
637 2095104 : return libMesh::RealGradient( value );
638 :
639 25141248 : case 1:
640 2095104 : return libMesh::RealGradient( Real(0), value );
641 :
642 25141248 : case 2:
643 2095104 : return libMesh::RealGradient( Real(0), Real(0), value );
644 :
645 0 : default:
646 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
647 : }
648 :
649 : //dummy
650 : return libMesh::RealGradient();
651 : }
652 :
653 31850496 : template <> RealGradient FE<3,HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
654 : const unsigned int i, const unsigned int j,
655 : const Point & p,
656 : const bool add_p_level)
657 : {
658 31850496 : const Real value = FE<3,HIERARCHIC>::shape_deriv(elem, order, i/3, j, p, add_p_level);
659 :
660 31850496 : switch( i%3 )
661 : {
662 10616832 : case 0:
663 884736 : return libMesh::RealGradient( value );
664 :
665 10616832 : case 1:
666 884736 : return libMesh::RealGradient( Real(0), value );
667 :
668 10616832 : case 2:
669 884736 : return libMesh::RealGradient( Real(0), Real(0), value );
670 :
671 0 : default:
672 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
673 : }
674 :
675 : //dummy
676 : return libMesh::RealGradient();
677 : }
678 :
679 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
680 :
681 0 : template <> RealGradient FE<3,HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
682 : const unsigned int i, const unsigned int j,
683 : const Point & p,
684 : const bool add_p_level)
685 : {
686 0 : const Real value = FE<3,HIERARCHIC>::shape_second_deriv(elem, order, i/3, j, p, add_p_level);
687 :
688 0 : switch( i%3 )
689 : {
690 0 : case 0:
691 0 : return libMesh::RealGradient( value );
692 :
693 0 : case 1:
694 0 : return libMesh::RealGradient( Real(0), value );
695 :
696 0 : case 2:
697 0 : return libMesh::RealGradient( Real(0), Real(0), value );
698 :
699 0 : default:
700 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
701 : }
702 :
703 : //dummy
704 : return libMesh::RealGradient();
705 : }
706 :
707 : #endif
708 :
709 75423744 : template <> RealGradient FE<3,L2_HIERARCHIC_VEC>::shape(const Elem * elem, const Order order,
710 : const unsigned int i, const Point & p,
711 : const bool add_p_level)
712 : {
713 75423744 : return FE<3,HIERARCHIC_VEC>::shape(elem, order, i, p, add_p_level);
714 : }
715 :
716 31850496 : template <> RealGradient FE<3,L2_HIERARCHIC_VEC>::shape_deriv(const Elem * elem, const Order order,
717 : const unsigned int i, const unsigned int j,
718 : const Point & p,
719 : const bool add_p_level)
720 : {
721 31850496 : return FE<3,HIERARCHIC_VEC>::shape_deriv(elem, order, i, j, p, add_p_level);
722 : }
723 :
724 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
725 :
726 0 : template <> RealGradient FE<3,L2_HIERARCHIC_VEC>::shape_second_deriv(const Elem * elem, const Order order,
727 : const unsigned int i, const unsigned int j,
728 : const Point & p,
729 : const bool add_p_level)
730 : {
731 0 : return FE<3,HIERARCHIC_VEC>::shape_second_deriv(elem, order, i, j, p, add_p_level);
732 : }
733 :
734 : #endif
735 :
736 :
737 : // Instantiate n_dofs*() functions for every dimension
738 695780 : LIBMESH_DEFAULT_VEC_NDOFS(HIERARCHIC)
739 :
740 : // L2 Hierarchic has a simpler but non-default implementation, with
741 : // all their dofs on the element
742 0 : template <> unsigned int FE<0,L2_HIERARCHIC_VEC>::n_dofs(const ElemType t, const Order o) { return FE<0,L2_HIERARCHIC>::n_dofs(t,o); }
743 0 : template <> unsigned int FE<1,L2_HIERARCHIC_VEC>::n_dofs(const ElemType t, const Order o) { return FE<1,L2_HIERARCHIC>::n_dofs(t,o); }
744 0 : template <> unsigned int FE<2,L2_HIERARCHIC_VEC>::n_dofs(const ElemType t, const Order o) { return 2*FE<2,L2_HIERARCHIC>::n_dofs(t,o); }
745 0 : template <> unsigned int FE<3,L2_HIERARCHIC_VEC>::n_dofs(const ElemType t, const Order o) { return 3*FE<3,L2_HIERARCHIC>::n_dofs(t,o); }
746 :
747 0 : template <> unsigned int FE<0,L2_HIERARCHIC_VEC>::n_dofs(const Elem * e, const Order o) { return FE<0,L2_HIERARCHIC>::n_dofs(e,o); }
748 0 : template <> unsigned int FE<1,L2_HIERARCHIC_VEC>::n_dofs(const Elem * e, const Order o) { return FE<1,L2_HIERARCHIC>::n_dofs(e,o); }
749 1727408 : template <> unsigned int FE<2,L2_HIERARCHIC_VEC>::n_dofs(const Elem * e, const Order o) { return 2*FE<2,L2_HIERARCHIC>::n_dofs(e,o); }
750 2014541 : template <> unsigned int FE<3,L2_HIERARCHIC_VEC>::n_dofs(const Elem * e, const Order o) { return 3*FE<3,L2_HIERARCHIC>::n_dofs(e,o); }
751 :
752 0 : template <> unsigned int FE<0,L2_HIERARCHIC_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
753 0 : template <> unsigned int FE<1,L2_HIERARCHIC_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
754 2706875 : template <> unsigned int FE<2,L2_HIERARCHIC_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
755 5173420 : template <> unsigned int FE<3,L2_HIERARCHIC_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
756 :
757 0 : template <> unsigned int FE<0,L2_HIERARCHIC_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
758 0 : template <> unsigned int FE<1,L2_HIERARCHIC_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
759 338143 : template <> unsigned int FE<2,L2_HIERARCHIC_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
760 771030 : template <> unsigned int FE<3,L2_HIERARCHIC_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
761 :
762 0 : template <> unsigned int FE<0,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<0,L2_HIERARCHIC_VEC>::n_dofs(t, o); }
763 0 : template <> unsigned int FE<1,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<1,L2_HIERARCHIC_VEC>::n_dofs(t, o); }
764 0 : template <> unsigned int FE<2,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<2,L2_HIERARCHIC_VEC>::n_dofs(t, o); }
765 0 : template <> unsigned int FE<3,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<3,L2_HIERARCHIC_VEC>::n_dofs(t, o); }
766 :
767 0 : template <> unsigned int FE<0,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<0,L2_HIERARCHIC_VEC>::n_dofs(&e, o); }
768 0 : template <> unsigned int FE<1,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<1,L2_HIERARCHIC_VEC>::n_dofs(&e, o); }
769 469656 : template <> unsigned int FE<2,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<2,L2_HIERARCHIC_VEC>::n_dofs(&e, o); }
770 412365 : template <> unsigned int FE<3,L2_HIERARCHIC_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<3,L2_HIERARCHIC_VEC>::n_dofs(&e, o); }
771 :
772 : // Hierarchic FEMs are always C^0 continuous
773 0 : template <> FEContinuity FE<0,HIERARCHIC_VEC>::get_continuity() const { return C_ZERO; }
774 0 : template <> FEContinuity FE<1,HIERARCHIC_VEC>::get_continuity() const { return C_ZERO; }
775 7860 : template <> FEContinuity FE<2,HIERARCHIC_VEC>::get_continuity() const { return C_ZERO; }
776 0 : template <> FEContinuity FE<3,HIERARCHIC_VEC>::get_continuity() const { return C_ZERO; }
777 :
778 : // L2 Hierarchic FEMs are always discontinuous
779 0 : template <> FEContinuity FE<0,L2_HIERARCHIC_VEC>::get_continuity() const { return DISCONTINUOUS; }
780 0 : template <> FEContinuity FE<1,L2_HIERARCHIC_VEC>::get_continuity() const { return DISCONTINUOUS; }
781 0 : template <> FEContinuity FE<2,L2_HIERARCHIC_VEC>::get_continuity() const { return DISCONTINUOUS; }
782 0 : template <> FEContinuity FE<3,L2_HIERARCHIC_VEC>::get_continuity() const { return DISCONTINUOUS; }
783 :
784 : // Hierarchic FEMs are hierarchic
785 0 : template <> bool FE<0,HIERARCHIC_VEC>::is_hierarchic() const { return true; }
786 0 : template <> bool FE<1,HIERARCHIC_VEC>::is_hierarchic() const { return true; }
787 0 : template <> bool FE<2,HIERARCHIC_VEC>::is_hierarchic() const { return true; }
788 0 : template <> bool FE<3,HIERARCHIC_VEC>::is_hierarchic() const { return true; }
789 0 : template <> bool FE<0,L2_HIERARCHIC_VEC>::is_hierarchic() const { return true; }
790 0 : template <> bool FE<1,L2_HIERARCHIC_VEC>::is_hierarchic() const { return true; }
791 0 : template <> bool FE<2,L2_HIERARCHIC_VEC>::is_hierarchic() const { return true; }
792 0 : template <> bool FE<3,L2_HIERARCHIC_VEC>::is_hierarchic() const { return true; }
793 :
794 : // Hierarchic FEM shapes need reinit
795 0 : template <> bool FE<0,HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
796 0 : template <> bool FE<1,HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
797 32958 : template <> bool FE<2,HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
798 0 : template <> bool FE<3,HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
799 0 : template <> bool FE<0,L2_HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
800 0 : template <> bool FE<1,L2_HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
801 1481176 : template <> bool FE<2,L2_HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
802 1791376 : template <> bool FE<3,L2_HIERARCHIC_VEC>::shapes_need_reinit() const { return true; }
803 :
804 : // Methods for computing Hierarchic constraints. Note: we pass the
805 : // dimension as the last argument to the anonymous helper function.
806 : // Also note: we only need instantiations of this function for
807 : // Dim==2 and 3.
808 : #ifdef LIBMESH_ENABLE_AMR
809 : template <>
810 7860 : void FE<2,HIERARCHIC_VEC>::compute_constraints (DofConstraints & constraints,
811 : DofMap & dof_map,
812 : const unsigned int variable_number,
813 : const Elem * elem)
814 : { //libmesh_not_implemented();
815 7860 : FEVectorBase::compute_proj_constraints(constraints, dof_map, variable_number, elem);
816 7860 : }
817 :
818 : template <>
819 0 : void FE<3,HIERARCHIC_VEC>::compute_constraints (DofConstraints & constraints,
820 : DofMap & dof_map,
821 : const unsigned int variable_number,
822 : const Elem * elem)
823 : { //libmesh_not_implemented();
824 0 : FEVectorBase::compute_proj_constraints(constraints, dof_map, variable_number, elem);
825 0 : }
826 :
827 : template <>
828 0 : void FE<2,L2_HIERARCHIC_VEC>::compute_constraints (DofConstraints & constraints,
829 : DofMap & dof_map,
830 : const unsigned int variable_number,
831 : const Elem * elem)
832 : { //libmesh_not_implemented();
833 0 : FEVectorBase::compute_proj_constraints(constraints, dof_map, variable_number, elem);
834 0 : }
835 :
836 : template <>
837 0 : void FE<3,L2_HIERARCHIC_VEC>::compute_constraints (DofConstraints & constraints,
838 : DofMap & dof_map,
839 : const unsigned int variable_number,
840 : const Elem * elem)
841 : { //libmesh_not_implemented();
842 0 : FEVectorBase::compute_proj_constraints(constraints, dof_map, variable_number, elem);
843 0 : }
844 : #endif // LIBMESH_ENABLE_AMR
845 :
846 : } // namespace libMesh
|