-
Notifications
You must be signed in to change notification settings - Fork 292
/
Copy pathfe.h
1758 lines (1545 loc) · 60.2 KB
/
fe.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
// The libMesh Finite Element Library.
// Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// Lesser General Public License for more details.
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
#ifndef LIBMESH_FE_H
#define LIBMESH_FE_H
// Local includes
#include "libmesh/fe_base.h"
#include "libmesh/int_range.h"
#include "libmesh/libmesh.h"
// C++ includes
#include <cstddef>
namespace libMesh
{
// forward declarations
class DofConstraints;
class DofMap;
class QGauss;
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
class InfFE;
#endif
/**
* Most finite element types in libMesh are scalar-valued
*/
template <FEFamily T>
struct FEOutputType
{
typedef Real type;
};
/**
* Specialize for non-scalar-valued elements
*/
template<>
struct FEOutputType<LAGRANGE_VEC>
{
typedef RealVectorValue type;
};
template<>
struct FEOutputType<L2_LAGRANGE_VEC>
{
typedef RealVectorValue type;
};
template<>
struct FEOutputType<HIERARCHIC_VEC>
{
typedef RealVectorValue type;
};
template<>
struct FEOutputType<L2_HIERARCHIC_VEC>
{
typedef RealVectorValue type;
};
template<>
struct FEOutputType<NEDELEC_ONE>
{
typedef RealVectorValue type;
};
template<>
struct FEOutputType<MONOMIAL_VEC>
{
typedef RealVectorValue type;
};
template<>
struct FEOutputType<RAVIART_THOMAS>
{
typedef RealVectorValue type;
};
template<>
struct FEOutputType<L2_RAVIART_THOMAS>
{
typedef RealVectorValue type;
};
/**
* A specific instantiation of the \p FEBase class. This
* class is templated, and specific template instantiations
* will result in different Finite Element families. Full specialization
* of the template for specific dimensions(\p Dim) and families
* (\p T) provide support for specific finite element types.
* The use of templates allows for compile-time optimization,
* however it requires that the specific finite element family
* and dimension is also known at compile time. If this is
* too restricting for your application you can use the
* \p FEBase::build() member to create abstract (but still optimized)
* finite elements.
*
* \author Benjamin S. Kirk
* \date 2002-2007
* \brief Template class which generates the different FE families and orders.
*/
template <unsigned int Dim, FEFamily T>
class FE : public FEGenericBase<typename FEOutputType<T>::type>
{
public:
/**
* Constructor.
*/
explicit
FE(const FEType & fet);
typedef typename
FEGenericBase<typename FEOutputType<T>::type>::OutputShape
OutputShape;
/**
* \returns The value of the \f$ i^{th} \f$ shape function at
* point \p p. This method allows you to specify the dimension,
* element type, and order directly. This allows the method to
* be static.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static OutputShape shape(const ElemType t,
const Order o,
const unsigned int i,
const Point & p);
/**
* \returns The value of the \f$ i^{th} \f$ shape function at
* point \p p. This method allows you to specify the dimension,
* element type, and order directly. This allows the method to
* be static.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static OutputShape shape(const Elem * elem,
const Order o,
const unsigned int i,
const Point & p,
const bool add_p_level = true);
/**
* \returns The value of the \f$ i^{th} \f$ shape function at
* point \p p. This method allows you to specify the dimension and
* element type directly. The order is given by the FEType.
* This allows the method to be static.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static OutputShape shape(const FEType fet,
const Elem * elem,
const unsigned int i,
const Point & p,
const bool add_p_level = true);
/**
* Fills \p v with the values of the \f$ i^{th} \f$
* shape function, evaluated at all points p. You must specify
* element order directly. \p v should already be the appropriate
* size.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static void shapes(const Elem * elem,
const Order o,
const unsigned int i,
const std::vector<Point> & p,
std::vector<OutputShape> & v,
const bool add_p_level = true);
/**
* Fills \p v[i][qp] with the values of the \f$ i^{th} \f$
* shape functions, evaluated at all points in p. You must specify
* element order directly. \p v should already be the appropriate
* size.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static void all_shapes(const Elem * elem,
const Order o,
const std::vector<Point> & p,
std::vector<std::vector<OutputShape> > & v,
const bool add_p_level = true);
/**
* \returns The \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
* shape function at point \p p. This method allows you to
* specify the dimension, element type, and order directly.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static OutputShape shape_deriv(const ElemType t,
const Order o,
const unsigned int i,
const unsigned int j,
const Point & p);
/**
* \returns The \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
* shape function. You must specify element type, and order directly.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static OutputShape shape_deriv(const Elem * elem,
const Order o,
const unsigned int i,
const unsigned int j,
const Point & p,
const bool add_p_level = true);
/**
* \returns The \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
* shape function. You must specify element type, and order (via
* FEType) directly.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static OutputShape shape_deriv(const FEType fet,
const Elem * elem,
const unsigned int i,
const unsigned int j,
const Point & p,
const bool add_p_level = true);
/**
* Fills \p v with the \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
* shape function, evaluated at all points p. You must specify
* element order directly. \p v should already be the appropriate
* size.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static void shape_derivs(const Elem * elem,
const Order o,
const unsigned int i,
const unsigned int j,
const std::vector<Point> & p,
std::vector<OutputShape> & v,
const bool add_p_level = true);
/**
* Fills \p comps with dphidxi (and in higher dimensions, eta/zeta)
* derivative component values for all shape functions, evaluated at
* all points in p. You must specify element order directly.
* Output component arrays in \p comps should already be the
* appropriate size.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static void all_shape_derivs(const Elem * elem,
const Order o,
const std::vector<Point> & p,
std::vector<std::vector<OutputShape>> * comps[3],
const bool add_p_level = true);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
/**
* \returns The second \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
* shape function at the point \p p.
*
* \note Cross-derivatives are indexed according to:
* j = 0 ==> d^2 phi / dxi^2
* j = 1 ==> d^2 phi / dxi deta
* j = 2 ==> d^2 phi / deta^2
* j = 3 ==> d^2 phi / dxi dzeta
* j = 4 ==> d^2 phi / deta dzeta
* j = 5 ==> d^2 phi / dzeta^2
*
* \note Computing second derivatives is not currently supported for
* all element types: \f$ C^1 \f$ (Clough, Hermite and Subdivision),
* Lagrange, Hierarchic, L2_Hierarchic, and Monomial are supported.
* All other element types return an error when asked for second
* derivatives.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static OutputShape shape_second_deriv(const ElemType t,
const Order o,
const unsigned int i,
const unsigned int j,
const Point & p);
/**
* \returns The second \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
* shape function at the point \p p.
*
* \note Cross-derivatives are indexed according to:
* j = 0 ==> d^2 phi / dxi^2
* j = 1 ==> d^2 phi / dxi deta
* j = 2 ==> d^2 phi / deta^2
* j = 3 ==> d^2 phi / dxi dzeta
* j = 4 ==> d^2 phi / deta dzeta
* j = 5 ==> d^2 phi / dzeta^2
*
* \note Computing second derivatives is not currently supported for
* all element types: \f$ C^1 \f$ (Clough, Hermite and Subdivision),
* Lagrange, Hierarchic, L2_Hierarchic, and Monomial are supported.
* All other element types return an error when asked for second
* derivatives.
*
* On a p-refined element, \p o should be the base order of the
* element if \p add_p_level is left \p true, or can be the base
* order of the element if \p add_p_level is set to \p false.
*/
static OutputShape shape_second_deriv(const Elem * elem,
const Order o,
const unsigned int i,
const unsigned int j,
const Point & p,
const bool add_p_level = true);
/**
* \returns The second \f$ j^{th} \f$ derivative of the \f$ i^{th} \f$
* shape function at the point \p p.
*
* \note Cross-derivatives are indexed according to:
* j = 0 ==> d^2 phi / dxi^2
* j = 1 ==> d^2 phi / dxi deta
* j = 2 ==> d^2 phi / deta^2
* j = 3 ==> d^2 phi / dxi dzeta
* j = 4 ==> d^2 phi / deta dzeta
* j = 5 ==> d^2 phi / dzeta^2
*
* \note Computing second derivatives is not currently supported for
* all element types: \f$ C^1 \f$ (Clough, Hermite and Subdivision),
* Lagrange, Hierarchic, L2_Hierarchic, and Monomial are supported.
* All other element types return an error when asked for second
* derivatives.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static OutputShape shape_second_deriv(const FEType fet,
const Elem * elem,
const unsigned int i,
const unsigned int j,
const Point & p,
const bool add_p_level = true);
#endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
/**
* Build the nodal soln from the element soln.
* This is the solution that will be plotted.
*
* On a p-refined element, \p o should be the base order of the element.
*/
static void nodal_soln(const Elem * elem, const Order o,
const std::vector<Number> & elem_soln,
std::vector<Number> & nodal_soln,
bool add_p_level = true,
const unsigned vdim = 1);
/**
* Build the nodal soln on one side from the (full) element soln.
* This is the solution that will be plotted on side-elements.
*
* On a p-refined element, \p o should be the base order of the element.
*/
static void side_nodal_soln(const Elem * elem, const Order o,
const unsigned int side,
const std::vector<Number> & elem_soln,
std::vector<Number> & nodal_soln_on_side,
bool add_p_level = true,
const unsigned vdim = 1);
/**
* \returns The number of shape functions associated with
* this finite element.
*/
virtual unsigned int n_shape_functions () const override;
/**
* \returns The number of shape functions associated with
* a finite element of type \p t and approximation order \p o.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method does not support all finite element types; e.g. for an
* arbitrary polygon or polyhedron type the number of shape
* functions may depend on an individual element and not just its
* type.
*/
static unsigned int n_shape_functions (const ElemType t,
const Order o)
{ return FE<Dim,T>::n_dofs (t,o); }
/**
* \returns The number of shape functions associated with this
* finite element.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method does not support all finite element types; e.g. for an
* arbitrary polygon or polyhedron type the number of shape
* functions may depend on an individual element and not just its
* type.
*/
static unsigned int n_dofs(const ElemType t,
const Order o);
/**
* \returns The number of shape functions associated with this
* finite element.
*
* On a p-refined element, \p o should be the total order of the element.
*
* \p e should only be a null pointer if using a FE family like
* SCALAR that has degrees of freedom independent of any element.
*/
static unsigned int n_dofs(const Elem * e,
const Order o);
/**
* \returns The number of dofs at node \p n for a finite element
* of type \p t and order \p o.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method does not support all finite element types; e.g. for an
* arbitrary polygon or polyhedron type the meaning of a node index
* \p n may depend on an individual element and not just its type.
*/
static unsigned int n_dofs_at_node(const ElemType t,
const Order o,
const unsigned int n);
/**
* \returns The number of dofs at node \p n for a finite element
* of type \p t and order \p o.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static unsigned int n_dofs_at_node(const Elem & e,
const Order o,
const unsigned int n);
/**
* \returns The number of dofs interior to the element,
* not associated with any interior nodes.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method may not support all finite element types, e.g. higher
* order polygons or polyhedra may differ from element to element.
*/
static unsigned int n_dofs_per_elem(const ElemType t,
const Order o);
/**
* \returns The number of dofs interior to the element,
* not associated with any interior nodes.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static unsigned int n_dofs_per_elem(const Elem & e,
const Order o);
/**
* \returns The continuity level of the finite element.
*/
virtual FEContinuity get_continuity() const override;
/**
* \returns \p true if the finite element's higher order shape functions are
* hierarchic
*/
virtual bool is_hierarchic() const override;
/**
* Fills the vector di with the local degree of freedom indices
* associated with side \p s of element \p elem
*
* On a p-refined element, \p o should be the base order of the element.
*/
static void dofs_on_side(const Elem * const elem,
const Order o,
unsigned int s,
std::vector<unsigned int> & di,
bool add_p_level=true);
/**
* Fills the vector di with the local degree of freedom indices
* associated with edge \p e of element \p elem
*
* On a p-refined element, \p o should be the base order of the element.
*/
static void dofs_on_edge(const Elem * const elem,
const Order o,
unsigned int e,
std::vector<unsigned int> & di,
bool add_p_level=true);
static Point inverse_map (const Elem * elem,
const Point & p,
const Real tolerance = TOLERANCE,
const bool secure = true)
{
// libmesh_deprecated(); // soon
return FEMap::inverse_map(Dim, elem, p, tolerance, secure, secure);
}
static void inverse_map (const Elem * elem,
const std::vector<Point> & physical_points,
std::vector<Point> & reference_points,
const Real tolerance = TOLERANCE,
const bool secure = true)
{
// libmesh_deprecated(); // soon
FEMap::inverse_map(Dim, elem, physical_points, reference_points,
tolerance, secure, secure);
}
/**
* This is at the core of this class. Use this for each
* new element in the mesh. Reinitializes all the physical
* element-dependent data based on the current element
* \p elem. By default the shape functions and associated
* data are computed at the quadrature points specified
* by the quadrature rule \p qrule, but may be any points
* specified on the reference element specified in the optional
* argument \p pts.
*/
virtual void reinit (const Elem * elem,
const std::vector<Point> * const pts = nullptr,
const std::vector<Real> * const weights = nullptr) override;
/**
* This re-computes the dual shape function coefficients.
* The dual shape coefficients are utilized when calculating dual shape functions.
*/
virtual void reinit_dual_shape_coeffs (const Elem * elem,
const std::vector<Point> & pts,
const std::vector<Real> & JxW) override;
/**
* This computes the default dual shape function coefficients.
* The dual shape coefficients are utilized when calculating dual shape functions.
*/
virtual void reinit_default_dual_shape_coeffs (const Elem * elem) override;
/**
* Reinitializes all the physical element-dependent data based on
* the \p side of \p face. The \p tolerance parameter is passed to
* the involved call to \p inverse_map(). By default the shape
* functions and associated data are computed at the quadrature
* points specified by the quadrature rule \p qrule, but may be any
* points specified on the reference \em side element specified in
* the optional argument \p pts.
*/
virtual void reinit (const Elem * elem,
const unsigned int side,
const Real tolerance = TOLERANCE,
const std::vector<Point> * const pts = nullptr,
const std::vector<Real> * const weights = nullptr) override;
/**
* Reinitializes all the physical element-dependent data based on
* the \p edge. The \p tolerance parameter is passed to the
* involved call to \p inverse_map(). By default the shape
* functions and associated data are computed at the quadrature
* points specified by the quadrature rule \p qrule, but may be any
* points specified on the reference \em side element specified in
* the optional argument \p pts.
*/
virtual void edge_reinit (const Elem * elem,
const unsigned int edge,
const Real tolerance = TOLERANCE,
const std::vector<Point> * const pts = nullptr,
const std::vector<Real> * const weights = nullptr) override;
/**
* Computes the reference space quadrature points on the side of
* an element based on the side quadrature points.
*/
virtual void side_map (const Elem * elem,
const Elem * side,
const unsigned int s,
const std::vector<Point> & reference_side_points,
std::vector<Point> & reference_points) override;
/**
* Computes the reference space quadrature points on the side of
* an element based on the edge quadrature points.
*/
virtual void edge_map (const Elem * elem,
const Elem * edge,
const unsigned int e,
const std::vector<Point> & reference_edge_points,
std::vector<Point> & reference_points);
/**
* Provides the class with the quadrature rule, which provides the
* locations (on a reference element) where the shape functions are
* to be calculated.
*/
virtual void attach_quadrature_rule (QBase * q) override;
/**
* \returns The total number of quadrature points. Call this
* to get an upper bound for the \p for loop in your simulation
* for matrix assembly of the current element.
*/
virtual unsigned int n_quadrature_points () const override;
#ifdef LIBMESH_ENABLE_AMR
/**
* Computes the constraint matrix contributions (for
* non-conforming adapted meshes) corresponding to
* variable number \p var_number, using element-specific
* optimizations if possible.
*/
static void compute_constraints (DofConstraints & constraints,
DofMap & dof_map,
const unsigned int variable_number,
const Elem * elem);
#endif // #ifdef LIBMESH_ENABLE_AMR
/**
* \returns \p true when the shape functions (for
* this \p FEFamily) depend on the particular
* element, and therefore needs to be re-initialized
* for each new element. \p false otherwise.
*/
virtual bool shapes_need_reinit() const override;
static Point map (const Elem * elem,
const Point & reference_point)
{
// libmesh_deprecated(); // soon
return FEMap::map(Dim, elem, reference_point);
}
static Point map_xi (const Elem * elem,
const Point & reference_point)
{
// libmesh_deprecated(); // soon
return FEMap::map_deriv(Dim, elem, 0, reference_point);
}
static Point map_eta (const Elem * elem,
const Point & reference_point)
{
// libmesh_deprecated(); // soon
return FEMap::map_deriv(Dim, elem, 1, reference_point);
}
static Point map_zeta (const Elem * elem,
const Point & reference_point)
{
// libmesh_deprecated(); // soon
return FEMap::map_deriv(Dim, elem, 2, reference_point);
}
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
/**
* make InfFE classes friends, so that these may access
* the private \p map, map_xyz methods
*/
template <unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
friend class InfFE;
#endif
protected:
/**
* Update the various member data fields \p phi,
* \p dphidxi, \p dphideta, \p dphidzeta, etc.
* for the current element. These data will be computed
* at the points \p qp, which are generally (but need not be)
* the quadrature points.
*/
virtual void init_shape_functions(const std::vector<Point> & qp,
const Elem * e);
/**
* A default implementation for all_shape_derivs
*/
static void default_all_shape_derivs (const Elem * elem,
const Order o,
const std::vector<Point> & p,
std::vector<std::vector<OutputShape>> * comps[3],
const bool add_p_level = true);
/**
* Init \p dual_phi and potentially \p dual_dphi, \p dual_d2phi
*/
void init_dual_shape_functions(unsigned int n_shapes, unsigned int n_qp);
#ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
/**
* Initialize the data fields for the base of an
* an infinite element.
*/
virtual void init_base_shape_functions(const std::vector<Point> & qp,
const Elem * e) override;
#endif
/**
* A default implementation for shapes
*/
static void default_shapes (const Elem * elem,
const Order o,
const unsigned int i,
const std::vector<Point> & p,
std::vector<OutputShape> & v,
const bool add_p_level = true)
{
libmesh_assert_equal_to(p.size(), v.size());
for (auto vi : index_range(v))
v[vi] = FE<Dim,T>::shape (elem, o, i, p[vi], add_p_level);
}
/**
* A default implementation for all_shapes
*/
static void default_all_shapes (const Elem * elem,
const Order o,
const std::vector<Point> & p,
std::vector<std::vector<OutputShape>> & v,
const bool add_p_level = true)
{
for (auto i : index_range(v))
{
libmesh_assert_equal_to ( p.size(), v[i].size() );
FE<Dim,T>::shapes (elem, o, i, p, v[i], add_p_level);
}
}
/**
* A default implementation for shape_derivs
*/
static void default_shape_derivs (const Elem * elem,
const Order o,
const unsigned int i,
const unsigned int j,
const std::vector<Point> & p,
std::vector<OutputShape> & v,
const bool add_p_level = true)
{
libmesh_assert_equal_to(p.size(), v.size());
for (auto vi : index_range(v))
v[vi] = FE<Dim,T>::shape_deriv (elem, o, i, j, p[vi], add_p_level);
}
/**
* A default implementation for side_nodal_soln
*/
static void default_side_nodal_soln(const Elem * elem, const Order o,
const unsigned int side,
const std::vector<Number> & elem_soln,
std::vector<Number> & nodal_soln_on_side,
bool add_p_level = true,
const unsigned vdim = 1);
/**
* An array of the node locations on the last
* element we computed on
*/
std::vector<Point> cached_nodes;
/**
* The last side and last edge we did a reinit on
*/
ElemType last_side;
ElemType last_edge;
};
/**
* Clough-Tocher finite elements. Still templated on the dimension,
* \p Dim.
*
* \author Roy Stogner
* \date 2004
*/
template <unsigned int Dim>
class FEClough : public FE<Dim,CLOUGH>
{
public:
/**
* Constructor. Creates a hierarchic finite element
* to be used in dimension \p Dim.
*/
explicit
FEClough(const FEType & fet) :
FE<Dim,CLOUGH> (fet)
{}
};
/**
* Hermite finite elements. Still templated on the dimension,
* \p Dim.
*
* \author Roy Stogner
* \date 2005
*/
template <unsigned int Dim>
class FEHermite : public FE<Dim,HERMITE>
{
public:
/**
* Constructor. Creates a hierarchic finite element
* to be used in dimension \p Dim.
*/
explicit
FEHermite(const FEType & fet) :
FE<Dim,HERMITE> (fet)
{}
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
/**
* 1D hermite functions on unit interval
*/
static Real hermite_raw_shape_second_deriv(const unsigned int basis_num,
const Real xi);
#endif
static Real hermite_raw_shape_deriv(const unsigned int basis_num,
const Real xi);
static Real hermite_raw_shape(const unsigned int basis_num,
const Real xi);
};
/**
* Subdivision finite elements.
*
* Template specialization prototypes are needed for calling from
* inside FESubdivision::init_shape_functions
*/
template <>
Real FE<2,SUBDIVISION>::shape(const Elem * elem,
const Order order,
const unsigned int i,
const Point & p,
const bool add_p_level);
template <>
Real FE<2,SUBDIVISION>::shape_deriv(const Elem * elem,
const Order order,
const unsigned int i,
const unsigned int j,
const Point & p,
const bool add_p_level);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
template <>
Real FE<2,SUBDIVISION>::shape_second_deriv(const Elem * elem,
const Order order,
const unsigned int i,
const unsigned int j,
const Point & p,
const bool add_p_level);
#endif
class FESubdivision : public FE<2,SUBDIVISION>
{
public:
/**
* Constructor. Creates a subdivision surface finite element.
* Currently only supported for two-dimensional meshes in
* three-dimensional space.
*/
FESubdivision(const FEType & fet);
/**
* This is at the core of this class. Use this for each new
* non-ghosted element in the mesh. Reinitializes all the physical
* element-dependent data based on the current element
* \p elem. By default the shape functions and associated
* data are computed at the quadrature points specified
* by the quadrature rule \p qrule, but may be any points
* specified on the reference element specified in the optional
* argument \p pts.
*/
virtual void reinit (const Elem * elem,
const std::vector<Point> * const pts = nullptr,
const std::vector<Real> * const weights = nullptr) override;
/**
* This prevents some compilers being confused by partially
* overriding this virtual function.
*/
virtual void reinit (const Elem *,
const unsigned int,
const Real = TOLERANCE,
const std::vector<Point> * const = nullptr,
const std::vector<Real> * const = nullptr) override
{ libmesh_not_implemented(); }
/**
* Provides the class with the quadrature rule, which provides the
* locations (on a reference element) where the shape functions are
* to be calculated.
*/
virtual void attach_quadrature_rule (QBase * q) override;
/**
* Update the various member data fields \p phi,
* \p dphidxi, \p dphideta, \p dphidzeta, etc.
* for the current element. These data will be computed
* at the points \p qp, which are generally (but need not be)
* the quadrature points.
*/
virtual void init_shape_functions(const std::vector<Point> & qp,
const Elem * elem) override;
/**
* \returns The value of the \f$ i^{th} \f$ of the 12 quartic
* box splines interpolating a regular Loop subdivision
* element, evaluated at the barycentric coordinates \p v,
* \p w.
*/
static Real regular_shape(const unsigned int i,
const Real v,
const Real w);
/**
* \returns The \f$ j^{th} \f$ derivative of the \f$ i^{th}
* \f$ of the 12 quartic box splines interpolating a regular
* Loop subdivision element, evaluated at the barycentric
* coordinates \p v, \p w.
*/
static Real regular_shape_deriv(const unsigned int i,
const unsigned int j,
const Real v,
const Real w);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
/**
* \returns The second \f$ j^{th} \f$ derivative of the
* \f$ i^{th} \f$ of the 12 quartic box splines interpolating
* a regular Loop subdivision element, evaluated at the
* barycentric coordinates \p v, \p w.
*/
static Real regular_shape_second_deriv(const unsigned int i,
const unsigned int j,
const Real v,
const Real w);
#endif // LIBMESH_ENABLE_SECOND_DERIVATIVE
/**
* Fills the vector \p weights with the weight coefficients
* of the Loop subdivision mask for evaluating the limit surface
* at a node explicitly. The size of \p weights will be
* 1 + \p valence, where \p valence is the number of neighbor
* nodes of the node where the limit surface is to be