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 : #ifndef LIBMESH_FE_INTERFACE_H
21 : #define LIBMESH_FE_INTERFACE_H
22 :
23 : // Local includes
24 : #include "libmesh/libmesh_common.h"
25 : #include "libmesh/vector_value.h"
26 :
27 : // C++ includes
28 : #include <map>
29 : #include <vector>
30 :
31 : namespace libMesh
32 : {
33 :
34 :
35 : // forward declarations
36 : class BoundaryInfo;
37 : class DofConstraints;
38 : class DofMap;
39 : class Elem;
40 : class FEType;
41 : class FEComputeData;
42 : class Point;
43 : class MeshBase;
44 : enum FEFamily : int;
45 : enum Order : int;
46 : enum FEFieldType : int;
47 : enum ElemType : int;
48 : enum FEContinuity : int;
49 :
50 : #ifdef LIBMESH_ENABLE_PERIODIC
51 : class PeriodicBoundaries;
52 : class PointLocatorBase;
53 : #endif
54 :
55 : /**
56 : * This class provides an encapsulated access to all static
57 : * public member functions of finite element classes.
58 : * Using this class, one need not worry about the correct
59 : * finite element class.
60 : *
61 : * \author Daniel Dreyer
62 : * \date 2002-2007
63 : * \brief Interface class which provides access to FE functions.
64 : */
65 : class FEInterface
66 : {
67 : private:
68 :
69 : /**
70 : * Empty constructor. Do not create an object of this type.
71 : */
72 : FEInterface();
73 :
74 : public:
75 :
76 : /**
77 : * Destructor.
78 : */
79 0 : virtual ~FEInterface() = default;
80 :
81 : #ifdef LIBMESH_ENABLE_DEPRECATED
82 : /**
83 : * \deprecated Call the version of this function taking an Elem* instead.
84 : */
85 : static unsigned int n_shape_functions(const unsigned int dim,
86 : const FEType & fe_t,
87 : const ElemType t);
88 : #endif // LIBMESH_ENABLE_DEPRECATED
89 :
90 : /**
91 : * \returns The number of shape functions associated with this
92 : * finite element \p elem of type \p fe_t.
93 : * Automatically decides which finite element class to use.
94 : *
95 : * On a p-refined element, \p fe_t.order should be the total order of the element.
96 : */
97 : static unsigned int n_shape_functions(const FEType & fe_t,
98 : const Elem * elem,
99 : const bool add_p_level = true);
100 :
101 : /**
102 : * Same as above, but ignores the elem->p_level() and uses the
103 : * specified extra_order instead.
104 : */
105 : static unsigned int n_shape_functions(const FEType & fe_t,
106 : const int extra_order,
107 : const Elem * elem);
108 :
109 : #ifdef LIBMESH_ENABLE_DEPRECATED
110 : /**
111 : * \deprecated Use n_dofs(const FEType &, Elem*) or n_dofs(const FEType &, int, Elem*) instead.
112 : */
113 : static unsigned int n_dofs(const unsigned int dim,
114 : const FEType & fe_t,
115 : const ElemType t);
116 :
117 : /**
118 : * Similar to the function above but takes an Elem * and accounts for
119 : * p-refinement internally, if any. This function is designed to prevent
120 : * users from needing to trick FEInterface::n_dofs() into giving them
121 : * the right number of dofs when working with p-refined elements. See,
122 : * e.g. FEInterface::compute_data().
123 : *
124 : * \deprecated Use n_dofs(const FEType &, Elem*) or n_dofs(const FEType &, int, Elem*) instead.
125 : */
126 : static unsigned int n_dofs(const unsigned int dim,
127 : const FEType & fe_t,
128 : const Elem * elem);
129 : #endif // LIBMESH_ENABLE_DEPRECATED
130 :
131 : /**
132 : * \returns The number of DOFs for \p elem for finite element type \p fe_t
133 : *
134 : * The p_level() of \p elem is accounted for internally by
135 : * increasing the Order of the passed-in FEType if \p add_p_level is true.
136 : */
137 : static unsigned int n_dofs(const FEType & fe_t,
138 : const Elem * elem,
139 : const bool add_p_level = true);
140 :
141 : /**
142 : * \returns The number of DOFs for \p elem for finite element type \p fe_t
143 : *
144 : * \note The p_level() of \elem is ignored and instead a total Order given
145 : * by fet_t.order + extra_order is used in determining the number of DOFs.
146 : */
147 : static unsigned int n_dofs(const FEType & fe_t,
148 : int extra_order,
149 : const Elem * elem);
150 :
151 : typedef unsigned int (*n_dofs_at_node_ptr) (const ElemType,
152 : const Order,
153 : const unsigned int);
154 :
155 : #ifdef LIBMESH_ENABLE_DEPRECATED
156 : /**
157 : * \returns The number of dofs at node n for a finite element
158 : * of type \p fe_t.
159 : * Automatically decides which finite element class to use.
160 : *
161 : * On a p-refined element, \p fe_t.order should be the total order of the element.
162 : *
163 : * \deprecated Call the version of n_dofs_at_node() taking an Elem *
164 : * instead, this one accounts for Elem::p_level() internally rather
165 : * than requiring the user to do it.
166 : */
167 : static unsigned int n_dofs_at_node(const unsigned int dim,
168 : const FEType & fe_t,
169 : const ElemType t,
170 : const unsigned int n);
171 :
172 : /**
173 : * \deprecated Use the version of this function that takes an Elem*
174 : * for consistency. The behavior is otherwise exactly the same,
175 : * since this function does not depend on the Elem::p_level().
176 : */
177 : static n_dofs_at_node_ptr
178 : n_dofs_at_node_function(const unsigned int dim,
179 : const FEType & fe_t);
180 : #endif // LIBMESH_ENABLE_DEPRECATED
181 :
182 : /**
183 : * \returns A function which evaluates n_dofs_at_node for the
184 : * requested FE type and element.
185 : */
186 : static n_dofs_at_node_ptr
187 : n_dofs_at_node_function(const FEType & fe_t,
188 : const Elem * elem);
189 :
190 : /**
191 : * \returns The number of dofs at node n for a finite element
192 : * of type \p fe_t. Accounts for Elem::p_level() internally if
193 : * \p add_p_level is true.
194 : */
195 : static unsigned int n_dofs_at_node(const FEType & fe_t,
196 : const Elem * elem,
197 : const unsigned int n,
198 : const bool add_p_level = true);
199 :
200 : /**
201 : * \returns The number of dofs at node n for a finite element
202 : * of type \p fe_t. Ignores Elem::p_level() and computes a total Order
203 : * given by fe_t.order + extra_order when determining the number of DOFs.
204 : */
205 : static unsigned int n_dofs_at_node(const FEType & fe_t,
206 : const int extra_order,
207 : const Elem * elem,
208 : const unsigned int n);
209 :
210 : #ifdef LIBMESH_ENABLE_DEPRECATED
211 : /**
212 : * \deprecated Call the version of this function that takes an Elem* instead.
213 : */
214 : static unsigned int n_dofs_per_elem(const unsigned int dim,
215 : const FEType & fe_t,
216 : const ElemType t);
217 : #endif // LIBMESH_ENABLE_DEPRECATED
218 :
219 : /**
220 : * \returns The number of dofs interior to the element,
221 : * not associated with any interior nodes.
222 : * Automatically decides which finite element class to use.
223 : *
224 : * On a p-refined element, \p fe_t.order should be the total order of the element.
225 : */
226 : static unsigned int n_dofs_per_elem(const FEType & fe_t,
227 : const Elem * elem,
228 : const bool add_p_level = true);
229 :
230 : /**
231 : * Same thing but internally elem->p_level() is ignored and extra_order is used instead.
232 : */
233 : static unsigned int n_dofs_per_elem(const FEType & fe_t,
234 : const int extra_order,
235 : const Elem * elem);
236 :
237 : /**
238 : * Fills the vector di with the local degree of freedom indices
239 : * associated with side \p s of element \p elem
240 : * Automatically decides which finite element class to use.
241 : *
242 : * On a p-refined element, \p fe_t.order should be the base order of the element.
243 : *
244 : * \todo For consistency with other FEInterface routines, this function
245 : * should be updated so that it does not take a \p dim argument.
246 : */
247 : static void dofs_on_side(const Elem * const elem,
248 : const unsigned int dim,
249 : const FEType & fe_t,
250 : unsigned int s,
251 : std::vector<unsigned int> & di,
252 : const bool add_p_level = true);
253 :
254 : /**
255 : * Fills the vector di with the local degree of freedom indices
256 : * associated with edge \p e of element \p elem
257 : * Automatically decides which finite element class to use.
258 : *
259 : * On a p-refined element, \p fe_t.order should be the base order of the element.
260 : *
261 : * \todo For consistency with other FEInterface routines, this function
262 : * should be updated so that it does not take a \p dim argument.
263 : */
264 : static void dofs_on_edge(const Elem * const elem,
265 : const unsigned int dim,
266 : const FEType & fe_t,
267 : unsigned int e,
268 : std::vector<unsigned int> & di,
269 : const bool add_p_level = true);
270 :
271 : /**
272 : * Build the nodal soln from the element soln.
273 : * This is the solution that will be plotted.
274 : * Automatically passes the request to the appropriate
275 : * finite element class member. To indicate that
276 : * results from this specific implementation of
277 : * \p nodal_soln should not be used, the vector
278 : * \p nodal_soln is returned empty.
279 : *
280 : * \note On a p-refined element, \p fe_t.order should be the base
281 : * order of the element. The Elem::p_level(), if any, is accounted
282 : * for internally by this routine.
283 : *
284 : * \todo For consistency with other FEInterface routines, this function
285 : * should be updated so that it does not take a \p dim argument.
286 : */
287 : static void nodal_soln(const unsigned int dim,
288 : const FEType & fe_t,
289 : const Elem * elem,
290 : const std::vector<Number> & elem_soln,
291 : std::vector<Number> & nodal_soln,
292 : const bool add_p_level = true,
293 : const unsigned int vdim = 1);
294 :
295 :
296 : /**
297 : * Build the nodal soln on one side from the (full) element soln.
298 : * This is the solution that will be plotted on side-elements.
299 : *
300 : * \note On a p-refined element, \p fe_t.order should be the base
301 : * order of the element. The Elem::p_level(), if any, is accounted
302 : * for internally by this routine.
303 : */
304 : static void side_nodal_soln(const FEType & fe_t,
305 : const Elem * elem,
306 : const unsigned int side,
307 : const std::vector<Number> & elem_soln,
308 : std::vector<Number> & nodal_soln,
309 : const bool add_p_level = true,
310 : const unsigned int vdim = 1);
311 :
312 : #ifdef LIBMESH_ENABLE_DEPRECATED
313 : /**
314 : * This is now deprecated; use FEMap::map instead.
315 : */
316 : static Point map(unsigned int dim,
317 : const FEType & fe_t,
318 : const Elem * elem,
319 : const Point & p);
320 :
321 : /**
322 : * This is now deprecated; use FEMap::inverse_map instead.
323 : */
324 : static Point inverse_map (const unsigned int dim,
325 : const FEType & fe_t,
326 : const Elem * elem,
327 : const Point & p,
328 : const Real tolerance = TOLERANCE,
329 : const bool secure = true);
330 :
331 : /**
332 : * This is now deprecated; use FEMap::inverse_map instead.
333 : */
334 : static void inverse_map (const unsigned int dim,
335 : const FEType & fe_t,
336 : const Elem * elem,
337 : const std::vector<Point> & physical_points,
338 : std::vector<Point> & reference_points,
339 : const Real tolerance = TOLERANCE,
340 : const bool secure = true);
341 :
342 : /**
343 : * \returns \p true if the point \p p is located on the reference element
344 : * for element type t, \p false otherwise.
345 : *
346 : * Since we are doing floating point comparisons, the parameter
347 : * \p eps can be specified to indicate a tolerance. For example,
348 : * \f$ \xi \le 1 \f$ becomes \f$ \xi \le 1 + \epsilon \f$.
349 : *
350 : * \deprecated This method overload does not support all finite
351 : * element types; e.g. the reference element for an arbitrary
352 : * polygon or polyhedron type may differ from element to element.
353 : * Use \p Elem::on_reference_element() instead.
354 : */
355 : static bool on_reference_element(const Point & p,
356 : const ElemType t,
357 : const Real eps=TOLERANCE);
358 :
359 : /**
360 : * \returns The value of the \f$ i^{th} \f$ shape function at
361 : * point \p p. This method allows you to specify the dimension,
362 : * element type, and order directly. Automatically passes the
363 : * request to the appropriate finite element class member.
364 : *
365 : * \note On a p-refined element, \p fe_t.order should be the total
366 : * order of the element.
367 : *
368 : * \deprecated Use the version of this function that accounts for
369 : * Elem::p_level() internally or the version which takes an
370 : * extra_order parameter.
371 : */
372 : static Real shape(const unsigned int dim,
373 : const FEType & fe_t,
374 : const ElemType t,
375 : const unsigned int i,
376 : const Point & p);
377 :
378 : /**
379 : * \returns The value of the \f$ i^{th} \f$ shape function at
380 : * point \p p. This method allows you to specify the dimension,
381 : * element type, and order directly. Automatically passes the
382 : * request to the appropriate finite element class member.
383 : *
384 : * \note On a p-refined element, \p fe_t.order should be the base
385 : * order of the element.
386 : *
387 : * \deprecated Use the version of this function that accounts for
388 : * Elem::p_level() internally or the version which takes an
389 : * extra_order parameter.
390 : */
391 : static Real shape(const unsigned int dim,
392 : const FEType & fe_t,
393 : const Elem * elem,
394 : const unsigned int i,
395 : const Point & p);
396 : #endif // LIBMESH_ENABLE_DEPRECATED
397 :
398 : /**
399 : * \returns The value of the \f$ i^{th} \f$ shape function at
400 : * point \p p.
401 : *
402 : * Non-deprecated version of the shape() function. The
403 : * Elem::p_level() is accounted for internally if \p add_p_level
404 : */
405 : static Real shape(const FEType & fe_t,
406 : const Elem * elem,
407 : const unsigned int i,
408 : const Point & p,
409 : const bool add_p_level = true);
410 :
411 : /**
412 : * \returns The value of the \f$ i^{th} \f$ shape function at
413 : * point \p p.
414 : *
415 : * Non-deprecated version of the shape() function. The
416 : * Elem::p_level() is ignored and extra_order is used instead.
417 : */
418 : static Real shape(const FEType & fe_t,
419 : int extra_order,
420 : const Elem * elem,
421 : const unsigned int i,
422 : const Point & p);
423 :
424 : #ifdef LIBMESH_ENABLE_DEPRECATED
425 : /**
426 : * \returns The value of the \f$ i^{th} \f$ shape function at
427 : * point \p p. This method allows you to specify the dimension,
428 : * element type, and order directly. Automatically passes the
429 : * request to the appropriate *scalar* finite element class member.
430 : *
431 : * \note On a p-refined element, \p fe_t.order should be the total
432 : * order of the element.
433 : *
434 : * \deprecated Use the version of this function that accounts for
435 : * Elem::p_level() internally or the version which takes an
436 : * extra_order parameter.
437 : */
438 : template<typename OutputType>
439 : static void shape(const unsigned int dim,
440 : const FEType & fe_t,
441 : const ElemType t,
442 : const unsigned int i,
443 : const Point & p,
444 : OutputType & phi);
445 :
446 : /**
447 : * \returns The value of the \f$ i^{th} \f$ shape function at
448 : * point \p p. This method allows you to specify the dimension,
449 : * element type, and order directly. Automatically passes the
450 : * request to the appropriate *scalar* finite element class member.
451 : *
452 : * \note On a p-refined element, \p fe_t.order should be the total
453 : * order of the element.
454 : *
455 : * \deprecated Use the version of this function that accounts for
456 : * Elem::p_level() internally or the version which takes an
457 : * extra_order parameter.
458 : */
459 : template<typename OutputType>
460 : static void shape(const unsigned int dim,
461 : const FEType & fe_t,
462 : const Elem * elem,
463 : const unsigned int i,
464 : const Point & p,
465 : OutputType & phi);
466 : #endif // LIBMESH_ENABLE_DEPRECATED
467 :
468 : /**
469 : * \returns The value of the \f$ i^{th} \f$ shape function at
470 : * point \p p. This method allows you to specify the dimension,
471 : * element type, and order directly. Automatically passes the
472 : * request to the appropriate *scalar* finite element class member.
473 : *
474 : * Non-deprecated version of templated shape() function that
475 : * accounts for Elem::p_level() internally.
476 : */
477 : template<typename OutputType>
478 : static void shape(const FEType & fe_t,
479 : const Elem * elem,
480 : const unsigned int i,
481 : const Point & p,
482 : OutputType & phi);
483 :
484 : /**
485 : * \returns The value of the \f$ i^{th} \f$ shape function at
486 : * point \p p. This method allows you to specify the dimension,
487 : * element type, and order directly. Automatically passes the
488 : * request to the appropriate *scalar* finite element class member.
489 : *
490 : * Non-deprecated version of templated shape() function that ignores
491 : * Elem::p_level() and instead uses extra_order internally.
492 : */
493 : template<typename OutputType>
494 : static void shape(const FEType & fe_t,
495 : int extra_order,
496 : const Elem * elem,
497 : const unsigned int i,
498 : const Point & p,
499 : OutputType & phi);
500 :
501 : /**
502 : * Fills \p phi with the values of the \f$ i^{th} \f$ shape function
503 : * at point \p p. This method allows you to specify the dimension,
504 : * element type, and order directly.
505 : *
506 : * \note Pass \p true for \p add_p_level if you want the Elem::p_level()
507 : * to be accounted for internally, pass false if you want fe_t.order to
508 : * be used instead.
509 : *
510 : * \todo To be consistent with the other non-deprecated FEInterface
511 : * routines, the shapes() and all_shapes() APIs should be updated so
512 : * that they do not take \p dim as a parameter. This is a relatively
513 : * large changeset with little benefit if we go the deprecation
514 : * route, so it would probably be cleaner to just break backwards
515 : * compatibility... these functions seem to mainly be used
516 : * internally by the library and changing them is unlikely to break
517 : * application codes.
518 : */
519 : template<typename OutputType>
520 : static void shapes(const unsigned int dim,
521 : const FEType & fe_t,
522 : const Elem * elem,
523 : const unsigned int i,
524 : const std::vector<Point> & p,
525 : std::vector<OutputType> & phi,
526 : const bool add_p_level = true);
527 :
528 : template<typename OutputType>
529 : static void all_shapes(const unsigned int dim,
530 : const FEType & fe_t,
531 : const Elem * elem,
532 : const std::vector<Point> & p,
533 : std::vector<std::vector<OutputType>> & phi,
534 : const bool add_p_level = true);
535 :
536 : /**
537 : * Typedef for pointer to a function that returns FE shape function values.
538 : * The \p p_level() of the passed-in \p elem is accounted for internally when
539 : * the \p add_p_level flag is set to true. For more information, see fe.h.
540 : */
541 : typedef Real (*shape_ptr) (const FEType fe_t,
542 : const Elem * elem,
543 : const unsigned int i,
544 : const Point & p,
545 : const bool add_p_level);
546 :
547 : #ifdef LIBMESH_ENABLE_DEPRECATED
548 : /**
549 : * \deprecated Call the version of this function taking an Elem* instead.
550 : */
551 : static shape_ptr
552 : shape_function(const unsigned int dim,
553 : const FEType & fe_t,
554 : const ElemType t);
555 : #endif // LIBMESH_ENABLE_DEPRECATED
556 :
557 : /**
558 : * \returns A function which evaluates shape for the
559 : * requested FE type and element.
560 : */
561 : static shape_ptr
562 : shape_function(const FEType & fe_t,
563 : const Elem * elem);
564 :
565 : #ifdef LIBMESH_ENABLE_DEPRECATED
566 : /**
567 : * \returns The \f$ j^{th} \f$ coordinate of the gradient of
568 : * the \f$ i^{th} \f$ shape function at point \p p.
569 : * This method allows you to specify the dimension,
570 : * element type, and order directly. Automatically passes the
571 : * request to the appropriate *scalar* finite element class member.
572 : *
573 : * \note On a p-refined element, \p fe_t.order should be the total
574 : * order of the element.
575 : *
576 : * \deprecated Call the version of this function taking an Elem* instead.
577 : */
578 : static Real shape_deriv(const unsigned int dim,
579 : const FEType & fe_t,
580 : const ElemType t,
581 : const unsigned int i,
582 : const unsigned int j,
583 : const Point & p);
584 :
585 : /**
586 : * \deprecated Call the version of this function taking an Elem* instead.
587 : */
588 : static Real shape_deriv (const unsigned int dim,
589 : const FEType & fe_t,
590 : const Elem *elem,
591 : const unsigned int i,
592 : const unsigned int j,
593 : const Point & p);
594 : #endif // LIBMESH_ENABLE_DEPRECATED
595 :
596 : /**
597 : * \returns The \f$ j^{th} \f$ coordinate of the gradient of
598 : * the \f$ i^{th} \f$ shape function at point \p p.
599 : * This method allows you to specify the dimension,
600 : * element, and order directly. Automatically passes the
601 : * request to the appropriate *scalar* finite element class member.
602 : *
603 : * \note On a p-refined element, \p fe_t.order should be the total
604 : * order of the element.
605 : */
606 : static Real shape_deriv(const FEType & fe_t,
607 : const Elem * elem,
608 : const unsigned int i,
609 : const unsigned int j,
610 : const Point & p);
611 :
612 : /**
613 : * Non-deprecated version of function above.
614 : */
615 : static Real shape_deriv(const FEType & fe_t,
616 : int extra_order,
617 : const Elem * elem,
618 : const unsigned int i,
619 : const unsigned int j,
620 : const Point & p);
621 :
622 : /**
623 : * Fills \p dphi with the derivatives of the \f$ i^{th} \f$ shape
624 : * function at point \p p in direction \p j.
625 : */
626 : template<typename OutputType>
627 : static void shape_derivs(const FEType & fe_t,
628 : const Elem * elem,
629 : const unsigned int i,
630 : const unsigned int j,
631 : const std::vector<Point> & p,
632 : std::vector<OutputType> & dphi,
633 : const bool add_p_level = true);
634 :
635 : template<typename OutputType>
636 : static void all_shape_derivs(const unsigned int dim,
637 : const FEType & fe_t,
638 : const Elem * elem,
639 : const std::vector<Point> & p,
640 : std::vector<std::vector<OutputType>> * comps[3],
641 : const bool add_p_level = true);
642 :
643 : /**
644 : * Typedef for pointer to a function that returns FE shape function derivative values.
645 : * The \p p_level() of the passed-in \p elem is accounted for internally when
646 : * the \p add_p_level flag is set to true. For more information, see fe.h.
647 : */
648 : typedef Real (*shape_deriv_ptr) (const FEType fet,
649 : const Elem * elem,
650 : const unsigned int i,
651 : const unsigned int j,
652 : const Point & p,
653 : const bool add_p_level);
654 :
655 : /**
656 : * \returns A function which evaluates shape for the
657 : * requested FE type and dimension.
658 : *
659 : * \deprecated Call the version of this function taking an Elem * instead.
660 : */
661 : static shape_deriv_ptr
662 : shape_deriv_function(const unsigned int dim,
663 : const FEType & fe_t,
664 : const ElemType t);
665 :
666 : /**
667 : * Non-deprecated version of the function above.
668 : */
669 : static shape_deriv_ptr
670 : shape_deriv_function(const FEType & fe_t,
671 : const Elem * elem);
672 :
673 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
674 :
675 : #ifdef LIBMESH_ENABLE_DEPRECATED
676 : /**
677 : * \returns The second \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
678 : * shape function at the point \p p.
679 : *
680 : * \note Cross-derivatives are indexed according to:
681 : * j = 0 ==> d^2 phi / dxi^2
682 : * j = 1 ==> d^2 phi / dxi deta
683 : * j = 2 ==> d^2 phi / deta^2
684 : * j = 3 ==> d^2 phi / dxi dzeta
685 : * j = 4 ==> d^2 phi / deta dzeta
686 : * j = 5 ==> d^2 phi / dzeta^2
687 : *
688 : * This method allows you to specify the dimension,
689 : * element type, and order directly. Automatically passes the
690 : * request to the appropriate *scalar* finite element class member.
691 : *
692 : * \note On a p-refined element, \p fe_t.order should be the total
693 : * order of the element.
694 : *
695 : * \deprecated Call version of this function which does not require
696 : * \p dim and takes an Elem * instead.
697 : */
698 : static Real shape_second_deriv(const unsigned int dim,
699 : const FEType & fe_t,
700 : const ElemType t,
701 : const unsigned int i,
702 : const unsigned int j,
703 : const Point & p);
704 :
705 : /**
706 : * \deprecated Call version of this function which does not require
707 : * \p dim and takes an Elem * instead.
708 : */
709 : static Real shape_second_deriv (const unsigned int dim,
710 : const FEType & fe_t,
711 : const Elem *elem,
712 : const unsigned int i,
713 : const unsigned int j,
714 : const Point & p);
715 : #endif // LIBMESH_ENABLE_DEPRECATED
716 :
717 : /**
718 : * \returns The second \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
719 : * shape function at the point \p p.
720 : *
721 : * \note Cross-derivatives are indexed according to:
722 : * j = 0 ==> d^2 phi / dxi^2
723 : * j = 1 ==> d^2 phi / dxi deta
724 : * j = 2 ==> d^2 phi / deta^2
725 : * j = 3 ==> d^2 phi / dxi dzeta
726 : * j = 4 ==> d^2 phi / deta dzeta
727 : * j = 5 ==> d^2 phi / dzeta^2
728 : *
729 : * \note On a p-refined element, \p fe_t.order should be the total
730 : * order of the element.
731 : */
732 : static Real shape_second_deriv(const FEType & fe_t,
733 : const Elem * elem,
734 : const unsigned int i,
735 : const unsigned int j,
736 : const Point & p);
737 :
738 : /**
739 : * Non-deprecated version of function above taking an \p extra_order parameter.
740 : */
741 : static Real shape_second_deriv(const FEType & fe_t,
742 : int extra_order,
743 : const Elem * elem,
744 : const unsigned int i,
745 : const unsigned int j,
746 : const Point & p);
747 :
748 : /**
749 : * Typedef for pointer to a function that returns FE shape function second derivative values.
750 : * The \p p_level() of the passed-in \p elem is accounted for internally when
751 : * the \p add_p_level flag is set to true. For more information, see fe.h.
752 : */
753 : typedef Real (*shape_second_deriv_ptr) (const FEType fet,
754 : const Elem * elem,
755 : const unsigned int i,
756 : const unsigned int j,
757 : const Point & p,
758 : const bool add_p_level);
759 :
760 : #ifdef LIBMESH_ENABLE_DEPRECATED
761 : /**
762 : * \deprecated Call the version of this function that takes an Elem * instead.
763 : */
764 : static shape_second_deriv_ptr
765 : shape_second_deriv_function(const unsigned int dim,
766 : const FEType & fe_t,
767 : const ElemType t);
768 : #endif // LIBMESH_ENABLE_DEPRECATED
769 :
770 : /**
771 : * \returns A function which evaluates shape for the
772 : * requested FE type and dimension.
773 : */
774 : static shape_second_deriv_ptr
775 : shape_second_deriv_function(const FEType & fe_t,
776 : const Elem * elem);
777 :
778 : #endif
779 :
780 : /**
781 : * Lets the appropriate child of \p FEBase compute the requested
782 : * data for the input specified in \p data, and sets the values in
783 : * \p data. See this as a generalization of \p shape(). With
784 : * infinite elements disabled, computes values for all shape
785 : * functions of \p elem evaluated at \p p.
786 : *
787 : * \note On a p-refined element, \p fe_t.order should be the base
788 : * order of the element.
789 : *
790 : * \todo For consistency with other FEInterface routines, this function
791 : * should be updated so that it does not take a \p dim argument.
792 : */
793 : static void compute_data(const unsigned int dim,
794 : const FEType & fe_t,
795 : const Elem * elem,
796 : FEComputeData & data);
797 :
798 : #ifdef LIBMESH_ENABLE_AMR
799 : /**
800 : * Computes the constraint matrix contributions (for
801 : * non-conforming adapted meshes) corresponding to
802 : * variable number \p var_number.
803 : */
804 : static void compute_constraints (DofConstraints & constraints,
805 : DofMap & dof_map,
806 : const unsigned int variable_number,
807 : const Elem * elem);
808 : #endif // #ifdef LIBMESH_ENABLE_AMR
809 :
810 : #ifdef LIBMESH_ENABLE_PERIODIC
811 : /**
812 : * Computes the constraint matrix contributions (for
813 : * periodic boundary conditions) corresponding to
814 : * variable number \p var_number.
815 : */
816 : static void compute_periodic_constraints (DofConstraints & constraints,
817 : DofMap & dof_map,
818 : const PeriodicBoundaries & boundaries,
819 : const MeshBase & mesh,
820 : const PointLocatorBase * point_locator,
821 : const unsigned int variable_number,
822 : const Elem * elem);
823 : #endif // #ifdef LIBMESH_ENABLE_PERIODIC
824 :
825 : /**
826 : * \returns The maximum polynomial degree that the given finite
827 : * element family can support on the given geometric element.
828 : */
829 : static unsigned int max_order (const FEType & fe_t,
830 : const ElemType & el_t);
831 :
832 : /**
833 : * \returns \p true if separate degrees of freedom must be allocated for
834 : * vertex DoFs and edge/face DoFs at a hanging node.
835 : */
836 : static bool extra_hanging_dofs (const FEType & fe_t);
837 :
838 : /**
839 : * \returns Whether the element is vector- or scalar-valued.
840 : */
841 : static FEFieldType field_type (const FEType & fe_type);
842 :
843 : /**
844 : * \returns Whether the element is vector- or scalar-valued.
845 : */
846 : static FEFieldType field_type (const FEFamily & fe_family);
847 :
848 : /**
849 : * \returns Whether the element's shape functions are orientation-dependent.
850 : */
851 : static bool orientation_dependent (const FEFamily & fe_family);
852 :
853 : /**
854 : * \returns The number of components of a vector-valued element.
855 : * Scalar-valued elements return 1.
856 : */
857 : static unsigned int n_vec_dim (const MeshBase & mesh,
858 : const FEType & fe_type);
859 :
860 : /**
861 : * Returns the input FEType's FEContinuity based on the underlying
862 : * FEFamily and potentially the Order, although we do not currently
863 : * support FEs with order-dependent continuity. These should exactly
864 : * match the FEBase::get_continuity() specializations/overrides for
865 : * the different FE types.
866 : */
867 : static FEContinuity get_continuity(const FEType & fe_type);
868 :
869 : /**
870 : * Returns whether or not the input FEType's higher-order shape
871 : * functions are always hierarchic
872 : */
873 : static bool is_hierarchic(const FEType & fe_type);
874 :
875 : private:
876 :
877 :
878 : /**
879 : * \returns \p true if \p et is an element to be processed by class
880 : * \p InfFE, false otherwise or when !LIBMESH_ENABLE_INFINITE_ELEMENTS.
881 : */
882 : static bool is_InfFE_elem(const ElemType et);
883 :
884 :
885 : #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
886 :
887 : /*
888 : * All these private members do the same as their public
889 : * counterparts, except for infinite elements. This disentangles
890 : * the calls to \p FE and \p InfFE.
891 : */
892 :
893 : #ifdef LIBMESH_ENABLE_DEPRECATED
894 : /**
895 : * \deprecated Call the version of ifem_n_shape_functions() which
896 : * takes a pointer-to-Elem instead.
897 : */
898 : static unsigned int ifem_n_shape_functions(const unsigned int dim,
899 : const FEType & fe_t,
900 : const ElemType t);
901 : #endif // LIBMESH_ENABLE_DEPRECATED
902 :
903 : static unsigned int ifem_n_shape_functions(const FEType & fe_t,
904 : const Elem * elem);
905 :
906 : #ifdef LIBMESH_ENABLE_DEPRECATED
907 : /**
908 : * \deprecated Call the version of ifem_n_dofs() which takes a
909 : * pointer-to-Elem instead.
910 : */
911 : static unsigned int ifem_n_dofs(const unsigned int dim,
912 : const FEType & fe_t,
913 : const ElemType t);
914 : #endif // LIBMESH_ENABLE_DEPRECATED
915 :
916 : static unsigned int ifem_n_dofs(const FEType & fe_t,
917 : const Elem * elem);
918 :
919 : #ifdef LIBMESH_ENABLE_DEPRECATED
920 : /**
921 : * \deprecated Call the version of ifem_n_dofs_at_node() which takes
922 : * a pointer-to-Elem instead.
923 : */
924 : static unsigned int ifem_n_dofs_at_node(const unsigned int dim,
925 : const FEType & fe_t,
926 : const ElemType t,
927 : const unsigned int n);
928 : #endif // LIBMESH_ENABLE_DEPRECATED
929 :
930 : static unsigned int ifem_n_dofs_at_node(const FEType & fe_t,
931 : const Elem * elem,
932 : const unsigned int n);
933 :
934 : #ifdef LIBMESH_ENABLE_DEPRECATED
935 : /**
936 : * \deprecated Call the version of ifem_n_dofs_per_elem() which
937 : * takes a pointer-to-Elem instead.
938 : */
939 : static unsigned int ifem_n_dofs_per_elem(const unsigned int dim,
940 : const FEType & fe_t,
941 : const ElemType t);
942 : #endif // LIBMESH_ENABLE_DEPRECATED
943 :
944 : static unsigned int ifem_n_dofs_per_elem(const FEType & fe_t,
945 : const Elem * elem);
946 :
947 : static void ifem_nodal_soln(const unsigned int dim,
948 : const FEType & fe_t,
949 : const Elem * elem,
950 : const std::vector<Number> & elem_soln,
951 : std::vector<Number> & nodal_soln);
952 :
953 : static Point ifem_map (const unsigned int dim,
954 : const FEType & fe_t,
955 : const Elem * elem,
956 : const Point & p);
957 :
958 : static Point ifem_inverse_map (const unsigned int dim,
959 : const FEType & fe_t,
960 : const Elem * elem,
961 : const Point & p,
962 : const Real tolerance = TOLERANCE,
963 : const bool secure = true);
964 :
965 : static void ifem_inverse_map (const unsigned int dim,
966 : const FEType & fe_t,
967 : const Elem * elem,
968 : const std::vector<Point> & physical_points,
969 : std::vector<Point> & reference_points,
970 : const Real tolerance = TOLERANCE,
971 : const bool secure = true);
972 :
973 : #ifdef LIBMESH_ENABLE_DEPRECATED
974 : /*
975 : * \deprecated This method overload may not support all finite
976 : * element types. Use \p Elem::on_reference_element() instead.
977 : */
978 : static bool ifem_on_reference_element(const Point & p,
979 : const ElemType t,
980 : const Real eps);
981 :
982 : /**
983 : * \deprecated Call version that takes a pointer-to-Elem and does
984 : * not require an explicit dim parameter instead.
985 : */
986 : static Real ifem_shape(const unsigned int dim,
987 : const FEType & fe_t,
988 : const ElemType t,
989 : const unsigned int i,
990 : const Point & p);
991 :
992 : /**
993 : * \deprecated Call version that takes a pointer-to-Elem and does
994 : * not require an explicit dim parameter instead.
995 : */
996 : static Real ifem_shape(const unsigned int dim,
997 : const FEType & fe_t,
998 : const Elem * elem,
999 : const unsigned int i,
1000 : const Point & p);
1001 : #endif // LIBMESH_ENABLE_DEPRECATED
1002 :
1003 : static Real ifem_shape(const FEType & fe_t,
1004 : const Elem * t,
1005 : const unsigned int i,
1006 : const Point & p);
1007 :
1008 : #ifdef LIBMESH_ENABLE_DEPRECATED
1009 : /**
1010 : * \deprecated Call version that takes a pointer-to-Elem and does
1011 : * not require an explicit dim parameter instead.
1012 : */
1013 : static Real ifem_shape_deriv(const unsigned int dim,
1014 : const FEType & fe_t,
1015 : const ElemType t,
1016 : const unsigned int i,
1017 : const unsigned int j,
1018 : const Point & p);
1019 :
1020 : /**
1021 : * \deprecated Call version that takes a pointer-to-Elem and does
1022 : * not require an explicit dim parameter instead.
1023 : */
1024 : static Real ifem_shape_deriv(const unsigned int dim,
1025 : const FEType & fe_t,
1026 : const Elem * elem,
1027 : const unsigned int i,
1028 : const unsigned int j,
1029 : const Point & p);
1030 : #endif // LIBMESH_ENABLE_DEPRECATED
1031 :
1032 : static Real ifem_shape_deriv(const FEType & fe_t,
1033 : const Elem * elem,
1034 : const unsigned int i,
1035 : const unsigned int j,
1036 : const Point & p);
1037 :
1038 : static void ifem_compute_data(const unsigned int dim,
1039 : const FEType & fe_t,
1040 : const Elem * elem,
1041 : FEComputeData & data);
1042 :
1043 : #endif
1044 :
1045 :
1046 : };
1047 :
1048 : } // namespace libMesh
1049 :
1050 : #endif // LIBMESH_FE_INTERFACE_H
|