MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_FE
13 #define MFEM_FE
14 
15 #include "../config/config.hpp"
16 #include "../general/array.hpp"
17 #include "../linalg/linalg.hpp"
18 #include "intrules.hpp"
19 #include "geom.hpp"
20 
21 namespace mfem
22 {
23 
24 // Base and derived classes for finite elements
25 
28 {
29 public:
30  enum
31  {
32  Pk, // polynomials of order k
33  Qk, // tensor products of polynomials of order k
34  rQk // refined tensor products of polynomials of order k
35  };
36 };
37 
38 class ElementTransformation;
39 class Coefficient;
40 class VectorCoefficient;
41 class KnotVector;
42 
45 {
46 protected:
49 
50 public:
52  enum { SCALAR, VECTOR };
53 
71  enum { VALUE, INTEGRAL, H_DIV, H_CURL };
72 
79  FiniteElement(int D, int G, int Do, int O, int F = FunctionSpace::Pk);
80 
82  int GetDim() const { return Dim; }
83 
85  int GetGeomType() const { return GeomType; }
86 
88  int GetDof() const { return Dof; }
89 
91  int GetOrder() const { return Order; }
92 
94  int Space() const { return FuncSpace; }
95 
96  int GetRangeType() const { return RangeType; }
97 
98  int GetMapType() const { return MapType; }
99 
103  virtual void CalcShape(const IntegrationPoint &ip,
104  Vector &shape) const = 0;
105 
110  virtual void CalcDShape(const IntegrationPoint &ip,
111  DenseMatrix &dshape) const = 0;
112 
113  const IntegrationRule & GetNodes() const { return Nodes; }
114 
115  // virtual functions for finite elements on vector spaces
116 
121  virtual void CalcVShape(const IntegrationPoint &ip,
122  DenseMatrix &shape) const;
123 
124  virtual void CalcVShape(ElementTransformation &Trans,
125  DenseMatrix &shape) const;
126 
130  virtual void CalcDivShape(const IntegrationPoint &ip,
131  Vector &divshape) const;
132 
137  virtual void CalcCurlShape(const IntegrationPoint &ip,
138  DenseMatrix &curl_shape) const;
139 
140  virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const;
141 
144  virtual void CalcHessian (const IntegrationPoint &ip,
145  DenseMatrix &h) const;
146 
150  virtual void GetLocalInterpolation (ElementTransformation &Trans,
151  DenseMatrix &I) const;
152 
156  virtual void Project (Coefficient &coeff,
157  ElementTransformation &Trans, Vector &dofs) const;
158 
162  virtual void Project (VectorCoefficient &vc,
163  ElementTransformation &Trans, Vector &dofs) const;
164 
167  virtual void ProjectDelta(int vertex, Vector &dofs) const;
168 
172  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
173  DenseMatrix &I) const;
174 
178  virtual void ProjectGrad(const FiniteElement &fe,
179  ElementTransformation &Trans,
180  DenseMatrix &grad) const;
181 
185  virtual void ProjectCurl(const FiniteElement &fe,
186  ElementTransformation &Trans,
187  DenseMatrix &curl) const;
188 
192  virtual void ProjectDiv(const FiniteElement &fe,
193  ElementTransformation &Trans,
194  DenseMatrix &div) const;
195 
196  virtual ~FiniteElement () { }
197 };
198 
200 {
201 protected:
203  DenseMatrix &I,
204  const NodalFiniteElement &fine_fe) const;
205 
206  void ProjectCurl_2D(const FiniteElement &fe,
207  ElementTransformation &Trans,
208  DenseMatrix &curl) const;
209 
210 #ifndef MFEM_THREAD_SAFE
211  mutable Vector c_shape;
212 #endif
213 
214 public:
215  NodalFiniteElement(int D, int G, int Do, int O,
216  int F = FunctionSpace::Pk) :
217 #ifdef MFEM_THREAD_SAFE
218  FiniteElement(D, G, Do, O, F)
219 #else
220  FiniteElement(D, G, Do, O, F), c_shape(Do)
221 #endif
222  { }
223 
224  void SetMapType(int M)
225  {
226  MFEM_VERIFY(M == VALUE || M == INTEGRAL, "unknown MapType");
227  MapType = M;
228  }
229 
231  DenseMatrix &I) const
232  { NodalLocalInterpolation (Trans, I, *this); }
233 
234  virtual void Project (Coefficient &coeff,
235  ElementTransformation &Trans, Vector &dofs) const;
236 
237  virtual void Project (VectorCoefficient &vc,
238  ElementTransformation &Trans, Vector &dofs) const;
239 
240  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
241  DenseMatrix &I) const;
242 
243  virtual void ProjectGrad(const FiniteElement &fe,
244  ElementTransformation &Trans,
245  DenseMatrix &grad) const;
246 
247  virtual void ProjectDiv(const FiniteElement &fe,
248  ElementTransformation &Trans,
249  DenseMatrix &div) const;
250 };
251 
252 
254 {
255 public:
256  PositiveFiniteElement(int D, int G, int Do, int O, int F = FunctionSpace::Pk)
257  : FiniteElement(D, G, Do, O, F) { }
259  virtual void Project(Coefficient &coeff,
260  ElementTransformation &Trans, Vector &dofs) const;
261  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
262  DenseMatrix &I) const;
263 };
264 
266 {
267  // Hide the scalar functions CalcShape and CalcDShape.
268 private:
270  virtual void CalcShape(const IntegrationPoint &ip,
271  Vector &shape) const;
272 
274  virtual void CalcDShape(const IntegrationPoint &ip,
275  DenseMatrix &dshape) const;
276 
277 protected:
278 #ifndef MFEM_THREAD_SAFE
279  mutable DenseMatrix J, Jinv;
281 #endif
282 
284  DenseMatrix &shape) const;
285 
287  DenseMatrix &shape) const;
288 
289  void Project_RT(const double *nk, const Array<int> &d2n,
291  Vector &dofs) const;
292 
293  void Project_RT(const double *nk, const Array<int> &d2n,
294  const FiniteElement &fe, ElementTransformation &Trans,
295  DenseMatrix &I) const;
296 
297  // rotated gradient in 2D
298  void ProjectGrad_RT(const double *nk, const Array<int> &d2n,
299  const FiniteElement &fe, ElementTransformation &Trans,
300  DenseMatrix &grad) const;
301 
302  void ProjectCurl_ND(const double *tk, const Array<int> &d2t,
303  const FiniteElement &fe, ElementTransformation &Trans,
304  DenseMatrix &curl) const;
305 
306  void ProjectCurl_RT(const double *nk, const Array<int> &d2n,
307  const FiniteElement &fe, ElementTransformation &Trans,
308  DenseMatrix &curl) const;
309 
310  void Project_ND(const double *tk, const Array<int> &d2t,
312  Vector &dofs) const;
313 
314  void Project_ND(const double *tk, const Array<int> &d2t,
315  const FiniteElement &fe, ElementTransformation &Trans,
316  DenseMatrix &I) const;
317 
318  void ProjectGrad_ND(const double *tk, const Array<int> &d2t,
319  const FiniteElement &fe, ElementTransformation &Trans,
320  DenseMatrix &grad) const;
321 
322  void LocalInterpolation_RT(const double *nk, const Array<int> &d2n,
323  ElementTransformation &Trans,
324  DenseMatrix &I) const;
325 
326  void LocalInterpolation_ND(const double *tk, const Array<int> &d2t,
327  ElementTransformation &Trans,
328  DenseMatrix &I) const;
329 
330 public:
331  VectorFiniteElement (int D, int G, int Do, int O, int M,
332  int F = FunctionSpace::Pk) :
333 #ifdef MFEM_THREAD_SAFE
334  FiniteElement(D, G, Do, O, F)
335  { RangeType = VECTOR; MapType = M; }
336 #else
337  FiniteElement(D, G, Do, O, F), Jinv(D), vshape(Do, D)
338  { RangeType = VECTOR; MapType = M; }
339 #endif
340 };
341 
343 {
344 public:
346 
347  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
348 
349  virtual void CalcDShape(const IntegrationPoint &ip,
350  DenseMatrix &dshape) const;
351 };
352 
355 {
356 public:
359 
363  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
364 
369  virtual void CalcDShape(const IntegrationPoint &ip,
370  DenseMatrix &dshape) const;
371 };
372 
375 {
376 public:
379 
383  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
384 
389  virtual void CalcDShape(const IntegrationPoint &ip,
390  DenseMatrix &dshape) const;
391  virtual void ProjectDelta(int vertex, Vector &dofs) const
392  { dofs = 0.0; dofs(vertex) = 1.0; }
393 };
394 
397 {
398 public:
401 
405  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
406 
411  virtual void CalcDShape(const IntegrationPoint &ip,
412  DenseMatrix &dshape) const;
413  virtual void CalcHessian (const IntegrationPoint &ip,
414  DenseMatrix &h) const;
415  virtual void ProjectDelta(int vertex, Vector &dofs) const
416  { dofs = 0.0; dofs(vertex) = 1.0; }
417  // { dofs = 1.0; }
418 };
419 
422 {
423 public:
425  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
426  virtual void CalcDShape(const IntegrationPoint &ip,
427  DenseMatrix &dshape) const;
428  virtual void ProjectDelta(int vertex, Vector &dofs) const;
429 };
430 
433 {
434 private:
435  static const double p[2];
436 
437 public:
439  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
440  virtual void CalcDShape(const IntegrationPoint &ip,
441  DenseMatrix &dshape) const;
442  virtual void ProjectDelta(int vertex, Vector &dofs) const;
443 };
444 
446 {
447 public:
449  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
450  virtual void CalcDShape(const IntegrationPoint &ip,
451  DenseMatrix &dshape) const;
452  virtual void ProjectDelta(int vertex, Vector &dofs) const
453  { dofs = 1.0; }
454 };
455 
458 {
459 public:
462 
466  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
467 
472  virtual void CalcDShape(const IntegrationPoint &ip,
473  DenseMatrix &dshape) const;
474 };
475 
477 {
478 public:
480  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
481  virtual void CalcDShape(const IntegrationPoint &ip,
482  DenseMatrix &dshape) const;
483 };
484 
487 {
488 public:
491 
495  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
496 
501  virtual void CalcDShape(const IntegrationPoint &ip,
502  DenseMatrix &dshape) const;
503 
504  virtual void CalcHessian (const IntegrationPoint &ip,
505  DenseMatrix &h) const;
506  virtual void ProjectDelta(int vertex, Vector &dofs) const;
507 };
508 
511 {
512 private:
513  static const double p[2];
514  DenseMatrix A;
515  mutable DenseMatrix D;
516  mutable Vector pol;
517 public:
519  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
520  virtual void CalcDShape(const IntegrationPoint &ip,
521  DenseMatrix &dshape) const;
522  // virtual void ProjectDelta(int vertex, Vector &dofs) const;
523 };
524 
527 {
528 public:
531 
535  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
536 
541  virtual void CalcDShape(const IntegrationPoint &ip,
542  DenseMatrix &dshape) const;
543  virtual void ProjectDelta(int vertex, Vector &dofs) const;
544 };
545 
547 {
548 public:
550  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
551  virtual void CalcDShape(const IntegrationPoint &ip,
552  DenseMatrix &dshape) const;
553  virtual void GetLocalInterpolation(ElementTransformation &Trans,
554  DenseMatrix &I) const;
556  virtual void Project(Coefficient &coeff, ElementTransformation &Trans,
557  Vector &dofs) const;
558  virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans,
559  Vector &dofs) const;
560  virtual void ProjectDelta(int vertex, Vector &dofs) const
561  { dofs = 0.; dofs(vertex) = 1.; }
562 };
563 
566 {
567 public:
569  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
570  virtual void CalcDShape(const IntegrationPoint &ip,
571  DenseMatrix &dshape) const;
572  // virtual void ProjectDelta(int vertex, Vector &dofs) const { dofs = 1.; }
573 };
574 
576 {
577 public:
579  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
580  virtual void CalcDShape(const IntegrationPoint &ip,
581  DenseMatrix &dshape) const;
582  virtual void CalcHessian (const IntegrationPoint &ip,
583  DenseMatrix &h) const;
584 };
585 
587 {
588 public:
590 
591  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
592 
593  virtual void CalcDShape(const IntegrationPoint &ip,
594  DenseMatrix &dshape) const;
595 };
596 
598 {
599 public:
601 
602  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
603 
604  virtual void CalcDShape(const IntegrationPoint &ip,
605  DenseMatrix &dshape) const;
606 
607  virtual void CalcHessian (const IntegrationPoint &ip,
608  DenseMatrix &h) const;
609 };
610 
613 {
614 public:
617 
618  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
619 
620  virtual void CalcDShape(const IntegrationPoint &ip,
621  DenseMatrix &dshape) const;
622 };
623 
626 {
627 public:
630 
632  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
633 
635  virtual void CalcDShape(const IntegrationPoint &ip,
636  DenseMatrix &dshape) const;
637  virtual void ProjectDelta(int vertex, Vector &dofs) const
638  { dofs(0) = 1.0; }
639 };
640 
641 
643 {
644 public:
646  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
647  virtual void CalcDShape(const IntegrationPoint &ip,
648  DenseMatrix &dshape) const;
649  virtual void ProjectDelta(int vertex, Vector &dofs) const
650  { dofs(0) = 1.0; }
651 };
652 
653 
656 {
657 public:
660 
664  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
665 
670  virtual void CalcDShape(const IntegrationPoint &ip,
671  DenseMatrix &dshape) const;
672 
673  virtual void ProjectDelta(int vertex, Vector &dofs) const
674  { dofs = 0.0; dofs(vertex) = 1.0; }
675 
676  virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const;
677 };
678 
681 {
682 public:
685 
686  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
687 
688  virtual void CalcDShape(const IntegrationPoint &ip,
689  DenseMatrix &dshape) const;
690 };
691 
694 {
695 public:
698 
702  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
703 
708  virtual void CalcDShape(const IntegrationPoint &ip,
709  DenseMatrix &dshape) const;
710 
711  virtual void ProjectDelta(int vertex, Vector &dofs) const
712  { dofs = 0.0; dofs(vertex) = 1.0; }
713 };
714 
717 {
718 public:
720  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
721  virtual void CalcDShape(const IntegrationPoint &ip,
722  DenseMatrix &dshape) const;
723  virtual void ProjectDelta(int vertex, Vector &dofs) const
724  { dofs = 1.0; }
725 };
726 
729 {
730 public:
732  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
733  virtual void CalcDShape(const IntegrationPoint &ip,
734  DenseMatrix &dshape) const;
735 };
736 
738 {
739 public:
740  P0SegmentFiniteElement(int Ord = 0);
741  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
742  virtual void CalcDShape(const IntegrationPoint &ip,
743  DenseMatrix &dshape) const;
744 };
745 
747 {
748 private:
749  static const double nk[3][2];
750 
751 public:
753 
754  virtual void CalcVShape(const IntegrationPoint &ip,
755  DenseMatrix &shape) const;
756 
757  virtual void CalcVShape(ElementTransformation &Trans,
758  DenseMatrix &shape) const
759  { CalcVShape_RT(Trans, shape); }
760 
761  virtual void CalcDivShape(const IntegrationPoint &ip,
762  Vector &divshape) const;
763 
764  virtual void GetLocalInterpolation (ElementTransformation &Trans,
765  DenseMatrix &I) const;
766 
768 
769  virtual void Project (VectorCoefficient &vc,
770  ElementTransformation &Trans, Vector &dofs) const;
771 };
772 
774 {
775 private:
776  static const double nk[4][2];
777 
778 public:
780 
781  virtual void CalcVShape(const IntegrationPoint &ip,
782  DenseMatrix &shape) const;
783 
784  virtual void CalcVShape(ElementTransformation &Trans,
785  DenseMatrix &shape) const
786  { CalcVShape_RT(Trans, shape); }
787 
788  virtual void CalcDivShape(const IntegrationPoint &ip,
789  Vector &divshape) const;
790 
791  virtual void GetLocalInterpolation (ElementTransformation &Trans,
792  DenseMatrix &I) const;
793 
795 
796  virtual void Project (VectorCoefficient &vc,
797  ElementTransformation &Trans, Vector &dofs) const;
798 };
799 
801 {
802 private:
803  static const double nk[8][2];
804 
805 public:
807 
808  virtual void CalcVShape(const IntegrationPoint &ip,
809  DenseMatrix &shape) const;
810 
811  virtual void CalcVShape(ElementTransformation &Trans,
812  DenseMatrix &shape) const
813  { CalcVShape_RT(Trans, shape); }
814 
815  virtual void CalcDivShape(const IntegrationPoint &ip,
816  Vector &divshape) const;
817 
818  virtual void GetLocalInterpolation (ElementTransformation &Trans,
819  DenseMatrix &I) const;
820 
822 
823  virtual void Project (VectorCoefficient &vc,
824  ElementTransformation &Trans, Vector &dofs) const;
825 };
826 
828 {
829 private:
830  static const double nk[12][2];
831 
832 public:
834 
835  virtual void CalcVShape(const IntegrationPoint &ip,
836  DenseMatrix &shape) const;
837 
838  virtual void CalcVShape(ElementTransformation &Trans,
839  DenseMatrix &shape) const
840  { CalcVShape_RT(Trans, shape); }
841 
842  virtual void CalcDivShape(const IntegrationPoint &ip,
843  Vector &divshape) const;
844 
845  virtual void GetLocalInterpolation (ElementTransformation &Trans,
846  DenseMatrix &I) const;
847 
849 
850  virtual void Project (VectorCoefficient &vc,
851  ElementTransformation &Trans, Vector &dofs) const;
852 };
853 
855 {
856 private:
857  static const double M[15][15];
858 public:
860 
861  virtual void CalcVShape(const IntegrationPoint &ip,
862  DenseMatrix &shape) const;
863 
864  virtual void CalcVShape(ElementTransformation &Trans,
865  DenseMatrix &shape) const
866  { CalcVShape_RT(Trans, shape); }
867 
868  virtual void CalcDivShape(const IntegrationPoint &ip,
869  Vector &divshape) const;
870 };
871 
873 {
874 private:
875  static const double nk[24][2];
876  static const double pt[4];
877  static const double dpt[3];
878 
879 public:
881 
882  virtual void CalcVShape(const IntegrationPoint &ip,
883  DenseMatrix &shape) const;
884 
885  virtual void CalcVShape(ElementTransformation &Trans,
886  DenseMatrix &shape) const
887  { CalcVShape_RT(Trans, shape); }
888 
889  virtual void CalcDivShape(const IntegrationPoint &ip,
890  Vector &divshape) const;
891 
892  virtual void GetLocalInterpolation (ElementTransformation &Trans,
893  DenseMatrix &I) const;
894 
896 
897  virtual void Project (VectorCoefficient &vc,
898  ElementTransformation &Trans, Vector &dofs) const;
899 };
900 
903 {
904 public:
906  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
907  virtual void CalcDShape(const IntegrationPoint &ip,
908  DenseMatrix &dshape) const;
909 };
910 
913 {
914 public:
916  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
917  virtual void CalcDShape(const IntegrationPoint &ip,
918  DenseMatrix &dshape) const;
919 };
920 
922 {
923 private:
924  Vector rwk;
925 #ifndef MFEM_THREAD_SAFE
926  mutable Vector rxxk;
927 #endif
928 public:
929  Lagrange1DFiniteElement (int degree);
930  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
931  virtual void CalcDShape(const IntegrationPoint &ip,
932  DenseMatrix &dshape) const;
933 };
934 
936 {
937 public:
939  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
940  virtual void CalcDShape(const IntegrationPoint &ip,
941  DenseMatrix &dshape) const;
942 };
943 
945 {
946 public:
948  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
949  virtual void CalcDShape(const IntegrationPoint &ip,
950  DenseMatrix &dshape) const;
951  virtual void ProjectDelta(int vertex, Vector &dofs) const
952  { dofs(0) = 1.0; }
953 };
954 
956 {
957 public:
959  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
960  virtual void CalcDShape(const IntegrationPoint &ip,
961  DenseMatrix &dshape) const;
962  virtual void ProjectDelta(int vertex, Vector &dofs) const
963  { dofs(0) = 1.0; }
964 };
965 
968 {
969 private:
971  int dof1d;
972  int *I, *J, *K;
973 #ifndef MFEM_THREAD_SAFE
974  mutable Vector shape1dx, shape1dy, shape1dz;
975  mutable DenseMatrix dshape1dx, dshape1dy, dshape1dz;
976 #endif
977 
978 public:
979  LagrangeHexFiniteElement (int degree);
980  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
981  virtual void CalcDShape(const IntegrationPoint &ip,
982  DenseMatrix &dshape) const;
984 };
985 
986 
989 {
990 public:
993 
997  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
998 
1003  virtual void CalcDShape(const IntegrationPoint &ip,
1004  DenseMatrix &dshape) const;
1005 };
1006 
1009 {
1010 public:
1013 
1017  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1018 
1023  virtual void CalcDShape(const IntegrationPoint &ip,
1024  DenseMatrix &dshape) const;
1025 };
1026 
1029 {
1030 public:
1033 
1034  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1035 
1036  virtual void CalcDShape(const IntegrationPoint &ip,
1037  DenseMatrix &dshape) const;
1038 };
1039 
1042 {
1043 public:
1046 
1050  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1051 
1056  virtual void CalcDShape(const IntegrationPoint &ip,
1057  DenseMatrix &dshape) const;
1058 };
1059 
1062 {
1063 public:
1066 
1070  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1071 
1076  virtual void CalcDShape(const IntegrationPoint &ip,
1077  DenseMatrix &dshape) const;
1078 };
1079 
1080 
1082 {
1083 private:
1084  static const double tk[12][3];
1085 
1086 public:
1088  virtual void CalcVShape(const IntegrationPoint &ip,
1089  DenseMatrix &shape) const;
1090  virtual void CalcVShape(ElementTransformation &Trans,
1091  DenseMatrix &shape) const
1092  { CalcVShape_ND(Trans, shape); }
1093  virtual void CalcCurlShape(const IntegrationPoint &ip,
1094  DenseMatrix &curl_shape) const;
1095  virtual void GetLocalInterpolation (ElementTransformation &Trans,
1096  DenseMatrix &I) const;
1097  using FiniteElement::Project;
1098  virtual void Project (VectorCoefficient &vc,
1099  ElementTransformation &Trans, Vector &dofs) const;
1100 };
1101 
1102 
1104 {
1105 private:
1106  static const double tk[6][3];
1107 
1108 public:
1110  virtual void CalcVShape(const IntegrationPoint &ip,
1111  DenseMatrix &shape) const;
1112  virtual void CalcVShape(ElementTransformation &Trans,
1113  DenseMatrix &shape) const
1114  { CalcVShape_ND(Trans, shape); }
1115  virtual void CalcCurlShape(const IntegrationPoint &ip,
1116  DenseMatrix &curl_shape) const;
1117  virtual void GetLocalInterpolation (ElementTransformation &Trans,
1118  DenseMatrix &I) const;
1119  using FiniteElement::Project;
1120  virtual void Project (VectorCoefficient &vc,
1121  ElementTransformation &Trans, Vector &dofs) const;
1122 };
1123 
1124 
1126 {
1127 private:
1128  static const double nk[6][3];
1129 
1130 public:
1132 
1133  virtual void CalcVShape(const IntegrationPoint &ip,
1134  DenseMatrix &shape) const;
1135 
1136  virtual void CalcVShape(ElementTransformation &Trans,
1137  DenseMatrix &shape) const
1138  { CalcVShape_RT(Trans, shape); }
1139 
1140  virtual void CalcDivShape(const IntegrationPoint &ip,
1141  Vector &divshape) const;
1142 
1143  virtual void GetLocalInterpolation (ElementTransformation &Trans,
1144  DenseMatrix &I) const;
1145 
1146  using FiniteElement::Project;
1147 
1148  virtual void Project (VectorCoefficient &vc,
1149  ElementTransformation &Trans, Vector &dofs) const;
1150 };
1151 
1152 
1154 {
1155 private:
1156  static const double nk[36][3];
1157 
1158 public:
1160 
1161  virtual void CalcVShape(const IntegrationPoint &ip,
1162  DenseMatrix &shape) const;
1163 
1164  virtual void CalcVShape(ElementTransformation &Trans,
1165  DenseMatrix &shape) const
1166  { CalcVShape_RT(Trans, shape); }
1167 
1168  virtual void CalcDivShape(const IntegrationPoint &ip,
1169  Vector &divshape) const;
1170 
1171  virtual void GetLocalInterpolation (ElementTransformation &Trans,
1172  DenseMatrix &I) const;
1173 
1174  using FiniteElement::Project;
1175 
1176  virtual void Project (VectorCoefficient &vc,
1177  ElementTransformation &Trans, Vector &dofs) const;
1178 };
1179 
1180 
1182 {
1183 private:
1184  static const double nk[4][3];
1185 
1186 public:
1188 
1189  virtual void CalcVShape(const IntegrationPoint &ip,
1190  DenseMatrix &shape) const;
1191 
1192  virtual void CalcVShape(ElementTransformation &Trans,
1193  DenseMatrix &shape) const
1194  { CalcVShape_RT(Trans, shape); }
1195 
1196  virtual void CalcDivShape(const IntegrationPoint &ip,
1197  Vector &divshape) const;
1198 
1199  virtual void GetLocalInterpolation (ElementTransformation &Trans,
1200  DenseMatrix &I) const;
1201 
1202  using FiniteElement::Project;
1203 
1204  virtual void Project (VectorCoefficient &vc,
1205  ElementTransformation &Trans, Vector &dofs) const;
1206 };
1207 
1208 
1210 {
1211 public:
1213  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1214  virtual void CalcDShape(const IntegrationPoint &ip,
1215  DenseMatrix &dshape) const;
1216 };
1217 
1218 
1219 class Poly_1D
1220 {
1221 public:
1222  class Basis
1223  {
1224  private:
1225  int mode; // 0 - use change of basis, O(p^2) Evals
1226  // 1 - use barycentric Lagrangian interpolation, O(p) Evals
1227  DenseMatrixInverse Ai;
1228  mutable Vector x, w;
1229 
1230  public:
1231  Basis(const int p, const double *nodes, const int _mode = 1);
1232  void Eval(const double x, Vector &u) const;
1233  void Eval(const double x, Vector &u, Vector &d) const;
1234  };
1235 
1236 private:
1237  Array<double *> open_pts, closed_pts;
1238  Array<Basis *> open_basis, closed_basis;
1239 
1240  static Array2D<int> binom;
1241 
1242  static void CalcMono(const int p, const double x, double *u);
1243  static void CalcMono(const int p, const double x, double *u, double *d);
1244 
1245  static void CalcLegendre(const int p, const double x, double *u);
1246  static void CalcLegendre(const int p, const double x, double *u, double *d);
1247 
1248  static void CalcChebyshev(const int p, const double x, double *u);
1249  static void CalcChebyshev(const int p, const double x, double *u, double *d);
1250 
1251 public:
1252  Poly_1D() { }
1253 
1254  static const int *Binom(const int p);
1255 
1256  const double *OpenPoints(const int p);
1257  const double *ClosedPoints(const int p);
1258 
1259  Basis &OpenBasis(const int p);
1260  Basis &ClosedBasis(const int p);
1261 
1262  // Evaluate the values of a hierarchical 1D basis at point x
1263  // hierarchical = k-th basis function is degree k polynomial
1264  static void CalcBasis(const int p, const double x, double *u)
1265  // { CalcMono(p, x, u); }
1266  // Bernstein basis is not hierarchical --> does not work for triangles
1267  // and tetrahedra
1268  // { CalcBernstein(p, x, u); }
1269  // { CalcLegendre(p, x, u); }
1270  { CalcChebyshev(p, x, u); }
1271 
1272  // Evaluate the values and derivatives of a hierarchical 1D basis at point x
1273  static void CalcBasis(const int p, const double x, double *u, double *d)
1274  // { CalcMono(p, x, u, d); }
1275  // { CalcBernstein(p, x, u, d); }
1276  // { CalcLegendre(p, x, u, d); }
1277  { CalcChebyshev(p, x, u, d); }
1278 
1279  // Evaluate a representation of a Delta function at point x
1280  static double CalcDelta(const int p, const double x)
1281  { return pow(x, (double) p); }
1282 
1283  static void UniformPoints(const int p, double *x);
1284  static void GaussPoints(const int p, double *x);
1285  static void GaussLobattoPoints(const int p, double *x);
1286  static void ChebyshevPoints(const int p, double *x);
1287 
1289  static void CalcBinomTerms(const int p, const double x, const double y,
1290  double *u);
1293  static void CalcBinomTerms(const int p, const double x, const double y,
1294  double *u, double *d);
1297  static void CalcDBinomTerms(const int p, const double x, const double y,
1298  double *d);
1299  static void CalcBernstein(const int p, const double x, double *u)
1300  { CalcBinomTerms(p, x, 1. - x, u); }
1301  static void CalcBernstein(const int p, const double x, double *u, double *d)
1302  { CalcBinomTerms(p, x, 1. - x, u, d); }
1303 
1304  ~Poly_1D();
1305 };
1306 
1307 extern Poly_1D poly1d;
1308 
1309 
1311 {
1312 private:
1313  Poly_1D::Basis &basis1d;
1314 #ifndef MFEM_THREAD_SAFE
1315  mutable Vector shape_x, dshape_x;
1316 #endif
1317 
1318 public:
1319  H1_SegmentElement(const int p);
1320  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1321  virtual void CalcDShape(const IntegrationPoint &ip,
1322  DenseMatrix &dshape) const;
1323  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1324 };
1325 
1326 
1328 {
1329 private:
1330  Poly_1D::Basis &basis1d;
1331 #ifndef MFEM_THREAD_SAFE
1332  mutable Vector shape_x, shape_y, dshape_x, dshape_y;
1333 #endif
1334  Array<int> dof_map;
1335 
1336 public:
1337  H1_QuadrilateralElement(const int p);
1338  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1339  virtual void CalcDShape(const IntegrationPoint &ip,
1340  DenseMatrix &dshape) const;
1341  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1342  const Array<int> &GetDofMap() const { return dof_map; }
1343 };
1344 
1345 
1347 {
1348 private:
1349  Poly_1D::Basis &basis1d;
1350 #ifndef MFEM_THREAD_SAFE
1351  mutable Vector shape_x, shape_y, shape_z, dshape_x, dshape_y, dshape_z;
1352 #endif
1353  Array<int> dof_map;
1354 
1355 public:
1356  H1_HexahedronElement(const int p);
1357  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1358  virtual void CalcDShape(const IntegrationPoint &ip,
1359  DenseMatrix &dshape) const;
1360  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1361  const Array<int> &GetDofMap() const { return dof_map; }
1362 };
1363 
1365 {
1366 private:
1367 #ifndef MFEM_THREAD_SAFE
1368  // This is to share scratch space between invocations, which helps
1369  // speed things up, but with OpenMP, we need one copy per thread.
1370  // Right now, we solve this by allocating this space within each function
1371  // call every time we call it. Alternatively, we should do some sort
1372  // thread private thing. Brunner, Jan 2014
1373  mutable Vector shape_x, dshape_x;
1374 #endif
1375 
1376 public:
1377  H1Pos_SegmentElement(const int p);
1378  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1379  virtual void CalcDShape(const IntegrationPoint &ip,
1380  DenseMatrix &dshape) const;
1381  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1382 };
1383 
1384 
1386 {
1387 private:
1388 #ifndef MFEM_THREAD_SAFE
1389  // See comment in H1Pos_SegmentElement
1390  mutable Vector shape_x, shape_y, dshape_x, dshape_y;
1391 #endif
1392  Array<int> dof_map;
1393 
1394 public:
1395  H1Pos_QuadrilateralElement(const int p);
1396  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1397  virtual void CalcDShape(const IntegrationPoint &ip,
1398  DenseMatrix &dshape) const;
1399  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1400 };
1401 
1402 
1404 {
1405 private:
1406 #ifndef MFEM_THREAD_SAFE
1407  // See comment in H1Pos_SegementElement.
1408  mutable Vector shape_x, shape_y, shape_z, dshape_x, dshape_y, dshape_z;
1409 #endif
1410  Array<int> dof_map;
1411 
1412 public:
1413  H1Pos_HexahedronElement(const int p);
1414  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1415  virtual void CalcDShape(const IntegrationPoint &ip,
1416  DenseMatrix &dshape) const;
1417  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1418 };
1419 
1420 
1422 {
1423 private:
1424 #ifndef MFEM_THREAD_SAFE
1425  mutable Vector shape_x, shape_y, shape_l, dshape_x, dshape_y, dshape_l, u;
1426  mutable DenseMatrix du;
1427 #endif
1428  DenseMatrixInverse Ti;
1429 
1430 public:
1431  H1_TriangleElement(const int p);
1432  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1433  virtual void CalcDShape(const IntegrationPoint &ip,
1434  DenseMatrix &dshape) const;
1435 };
1436 
1437 
1439 {
1440 private:
1441 #ifndef MFEM_THREAD_SAFE
1442  mutable Vector shape_x, shape_y, shape_z, shape_l;
1443  mutable Vector dshape_x, dshape_y, dshape_z, dshape_l, u;
1444  mutable DenseMatrix du;
1445 #endif
1446  DenseMatrixInverse Ti;
1447 
1448 public:
1449  H1_TetrahedronElement(const int p);
1450  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1451  virtual void CalcDShape(const IntegrationPoint &ip,
1452  DenseMatrix &dshape) const;
1453 };
1454 
1455 
1456 // TODO: define H1Pos_ FEs on triangle and tetrahedron
1457 
1458 
1460 {
1461 private:
1462  int type;
1463  Poly_1D::Basis *basis1d;
1464 #ifndef MFEM_THREAD_SAFE
1465  mutable Vector shape_x, dshape_x;
1466 #endif
1467 
1468 public:
1469  L2_SegmentElement(const int p, const int _type = 0);
1470  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1471  virtual void CalcDShape(const IntegrationPoint &ip,
1472  DenseMatrix &dshape) const;
1473  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1474 };
1475 
1476 
1478 {
1479 private:
1480 #ifndef MFEM_THREAD_SAFE
1481  mutable Vector shape_x, dshape_x;
1482 #endif
1483 
1484 public:
1485  L2Pos_SegmentElement(const int p);
1486  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1487  virtual void CalcDShape(const IntegrationPoint &ip,
1488  DenseMatrix &dshape) const;
1489  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1490 };
1491 
1492 
1494 {
1495 private:
1496  int type;
1497  Poly_1D::Basis *basis1d;
1498 #ifndef MFEM_THREAD_SAFE
1499  mutable Vector shape_x, shape_y, dshape_x, dshape_y;
1500 #endif
1501 
1502 public:
1503  L2_QuadrilateralElement(const int p, const int _type = 0);
1504  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1505  virtual void CalcDShape(const IntegrationPoint &ip,
1506  DenseMatrix &dshape) const;
1507  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1508  virtual void ProjectCurl(const FiniteElement &fe,
1509  ElementTransformation &Trans,
1510  DenseMatrix &curl) const
1511  { ProjectCurl_2D(fe, Trans, curl); }
1512 };
1513 
1514 
1516 {
1517 private:
1518 #ifndef MFEM_THREAD_SAFE
1519  mutable Vector shape_x, shape_y, dshape_x, dshape_y;
1520 #endif
1521 
1522 public:
1523  L2Pos_QuadrilateralElement(const int p);
1524  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1525  virtual void CalcDShape(const IntegrationPoint &ip,
1526  DenseMatrix &dshape) const;
1527  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1528 };
1529 
1530 
1532 {
1533 private:
1534  int type;
1535  Poly_1D::Basis *basis1d;
1536 #ifndef MFEM_THREAD_SAFE
1537  mutable Vector shape_x, shape_y, shape_z, dshape_x, dshape_y, dshape_z;
1538 #endif
1539 
1540 public:
1541  L2_HexahedronElement(const int p, const int _type = 0);
1542  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1543  virtual void CalcDShape(const IntegrationPoint &ip,
1544  DenseMatrix &dshape) const;
1545  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1546 };
1547 
1548 
1550 {
1551 private:
1552 #ifndef MFEM_THREAD_SAFE
1553  mutable Vector shape_x, shape_y, shape_z, dshape_x, dshape_y, dshape_z;
1554 #endif
1555 
1556 public:
1557  L2Pos_HexahedronElement(const int p);
1558  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1559  virtual void CalcDShape(const IntegrationPoint &ip,
1560  DenseMatrix &dshape) const;
1561  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1562 };
1563 
1564 
1566 {
1567 private:
1568  int type;
1569 #ifndef MFEM_THREAD_SAFE
1570  mutable Vector shape_x, shape_y, shape_l, dshape_x, dshape_y, dshape_l, u;
1571  mutable DenseMatrix du;
1572 #endif
1573  DenseMatrixInverse Ti;
1574 
1575 public:
1576  L2_TriangleElement(const int p, const int _type = 0);
1577  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1578  virtual void CalcDShape(const IntegrationPoint &ip,
1579  DenseMatrix &dshape) const;
1580  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1581  virtual void ProjectCurl(const FiniteElement &fe,
1582  ElementTransformation &Trans,
1583  DenseMatrix &curl) const
1584  { ProjectCurl_2D(fe, Trans, curl); }
1585 };
1586 
1587 
1589 {
1590 private:
1591 #ifndef MFEM_THREAD_SAFE
1592  mutable Vector dshape_1d;
1593 #endif
1594 
1595 public:
1596  L2Pos_TriangleElement(const int p);
1597  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1598  virtual void CalcDShape(const IntegrationPoint &ip,
1599  DenseMatrix &dshape) const;
1600  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1601 };
1602 
1603 
1605 {
1606 private:
1607  int type;
1608 #ifndef MFEM_THREAD_SAFE
1609  mutable Vector shape_x, shape_y, shape_z, shape_l;
1610  mutable Vector dshape_x, dshape_y, dshape_z, dshape_l, u;
1611  mutable DenseMatrix du;
1612 #endif
1613  DenseMatrixInverse Ti;
1614 
1615 public:
1616  L2_TetrahedronElement(const int p, const int _type = 0);
1617  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1618  virtual void CalcDShape(const IntegrationPoint &ip,
1619  DenseMatrix &dshape) const;
1620  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1621 };
1622 
1623 
1625 {
1626 private:
1627 #ifndef MFEM_THREAD_SAFE
1628  mutable Vector dshape_1d;
1629 #endif
1630 
1631 public:
1632  L2Pos_TetrahedronElement(const int p);
1633  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
1634  virtual void CalcDShape(const IntegrationPoint &ip,
1635  DenseMatrix &dshape) const;
1636  virtual void ProjectDelta(int vertex, Vector &dofs) const;
1637 };
1638 
1639 
1641 {
1642 private:
1643  static const double nk[8];
1644 
1645  Poly_1D::Basis &cbasis1d, &obasis1d;
1646 #ifndef MFEM_THREAD_SAFE
1647  mutable Vector shape_cx, shape_ox, shape_cy, shape_oy;
1648  mutable Vector dshape_cx, dshape_cy;
1649 #endif
1650  Array<int> dof_map, dof2nk;
1651 
1652 public:
1653  RT_QuadrilateralElement(const int p);
1654  virtual void CalcVShape(const IntegrationPoint &ip,
1655  DenseMatrix &shape) const;
1656  virtual void CalcVShape(ElementTransformation &Trans,
1657  DenseMatrix &shape) const
1658  { CalcVShape_RT(Trans, shape); }
1659  virtual void CalcDivShape(const IntegrationPoint &ip,
1660  Vector &divshape) const;
1662  DenseMatrix &I) const
1663  { LocalInterpolation_RT(nk, dof2nk, Trans, I); }
1664  using FiniteElement::Project;
1665  virtual void Project(VectorCoefficient &vc,
1666  ElementTransformation &Trans, Vector &dofs) const
1667  { Project_RT(nk, dof2nk, vc, Trans, dofs); }
1668  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
1669  DenseMatrix &I) const
1670  { Project_RT(nk, dof2nk, fe, Trans, I); }
1671  virtual void ProjectGrad(const FiniteElement &fe,
1672  ElementTransformation &Trans,
1673  DenseMatrix &grad) const
1674  { ProjectGrad_RT(nk, dof2nk, fe, Trans, grad); }
1675 };
1676 
1677 
1679 {
1680  static const double nk[18];
1681 
1682  Poly_1D::Basis &cbasis1d, &obasis1d;
1683 #ifndef MFEM_THREAD_SAFE
1684  mutable Vector shape_cx, shape_ox, shape_cy, shape_oy, shape_cz, shape_oz;
1685  mutable Vector dshape_cx, dshape_cy, dshape_cz;
1686 #endif
1687  Array<int> dof_map, dof2nk;
1688 
1689 public:
1690  RT_HexahedronElement(const int p);
1691  virtual void CalcVShape(const IntegrationPoint &ip,
1692  DenseMatrix &shape) const;
1693  virtual void CalcVShape(ElementTransformation &Trans,
1694  DenseMatrix &shape) const
1695  { CalcVShape_RT(Trans, shape); }
1696  virtual void CalcDivShape(const IntegrationPoint &ip,
1697  Vector &divshape) const;
1699  DenseMatrix &I) const
1700  { LocalInterpolation_RT(nk, dof2nk, Trans, I); }
1701  using FiniteElement::Project;
1702  virtual void Project(VectorCoefficient &vc,
1703  ElementTransformation &Trans, Vector &dofs) const
1704  { Project_RT(nk, dof2nk, vc, Trans, dofs); }
1705  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
1706  DenseMatrix &I) const
1707  { Project_RT(nk, dof2nk, fe, Trans, I); }
1708  virtual void ProjectCurl(const FiniteElement &fe,
1709  ElementTransformation &Trans,
1710  DenseMatrix &curl) const
1711  { ProjectCurl_RT(nk, dof2nk, fe, Trans, curl); }
1712 };
1713 
1714 
1716 {
1717  static const double nk[6], c;
1718 
1719 #ifndef MFEM_THREAD_SAFE
1720  mutable Vector shape_x, shape_y, shape_l;
1721  mutable Vector dshape_x, dshape_y, dshape_l;
1722  mutable DenseMatrix u;
1723  mutable Vector divu;
1724 #endif
1725  Array<int> dof2nk;
1726  DenseMatrixInverse Ti;
1727 
1728 public:
1729  RT_TriangleElement(const int p);
1730  virtual void CalcVShape(const IntegrationPoint &ip,
1731  DenseMatrix &shape) const;
1732  virtual void CalcVShape(ElementTransformation &Trans,
1733  DenseMatrix &shape) const
1734  { CalcVShape_RT(Trans, shape); }
1735  virtual void CalcDivShape(const IntegrationPoint &ip,
1736  Vector &divshape) const;
1738  DenseMatrix &I) const
1739  { LocalInterpolation_RT(nk, dof2nk, Trans, I); }
1740  using FiniteElement::Project;
1741  virtual void Project(VectorCoefficient &vc,
1742  ElementTransformation &Trans, Vector &dofs) const
1743  { Project_RT(nk, dof2nk, vc, Trans, dofs); }
1744  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
1745  DenseMatrix &I) const
1746  { Project_RT(nk, dof2nk, fe, Trans, I); }
1747  virtual void ProjectGrad(const FiniteElement &fe,
1748  ElementTransformation &Trans,
1749  DenseMatrix &grad) const
1750  { ProjectGrad_RT(nk, dof2nk, fe, Trans, grad); }
1751 };
1752 
1753 
1755 {
1756  static const double nk[12], c;
1757 
1758 #ifndef MFEM_THREAD_SAFE
1759  mutable Vector shape_x, shape_y, shape_z, shape_l;
1760  mutable Vector dshape_x, dshape_y, dshape_z, dshape_l;
1761  mutable DenseMatrix u;
1762  mutable Vector divu;
1763 #endif
1764  Array<int> dof2nk;
1765  DenseMatrixInverse Ti;
1766 
1767 public:
1768  RT_TetrahedronElement(const int p);
1769  virtual void CalcVShape(const IntegrationPoint &ip,
1770  DenseMatrix &shape) const;
1771  virtual void CalcVShape(ElementTransformation &Trans,
1772  DenseMatrix &shape) const
1773  { CalcVShape_RT(Trans, shape); }
1774  virtual void CalcDivShape(const IntegrationPoint &ip,
1775  Vector &divshape) const;
1777  DenseMatrix &I) const
1778  { LocalInterpolation_RT(nk, dof2nk, Trans, I); }
1779  using FiniteElement::Project;
1780  virtual void Project(VectorCoefficient &vc,
1781  ElementTransformation &Trans, Vector &dofs) const
1782  { Project_RT(nk, dof2nk, vc, Trans, dofs); }
1783  virtual void Project(const FiniteElement &fe, ElementTransformation &Trans,
1784  DenseMatrix &I) const
1785  { Project_RT(nk, dof2nk, fe, Trans, I); }
1786  virtual void ProjectCurl(const FiniteElement &fe,
1787  ElementTransformation &Trans,
1788  DenseMatrix &curl) const
1789  { ProjectCurl_RT(nk, dof2nk, fe, Trans, curl); }
1790 };
1791 
1792 
1794 {
1795  static const double tk[18];
1796 
1797  Poly_1D::Basis &cbasis1d, &obasis1d;
1798 #ifndef MFEM_THREAD_SAFE
1799  mutable Vector shape_cx, shape_ox, shape_cy, shape_oy, shape_cz, shape_oz;
1800  mutable Vector dshape_cx, dshape_cy, dshape_cz;
1801 #endif
1802  Array<int> dof_map, dof2tk;
1803 
1804 public:
1805  ND_HexahedronElement(const int p);
1806 
1807  virtual void CalcVShape(const IntegrationPoint &ip,
1808  DenseMatrix &shape) const;
1809 
1810  virtual void CalcVShape(ElementTransformation &Trans,
1811  DenseMatrix &shape) const
1812  { CalcVShape_ND(Trans, shape); }
1813 
1814  virtual void CalcCurlShape(const IntegrationPoint &ip,
1815  DenseMatrix &curl_shape) const;
1816 
1818  DenseMatrix &I) const
1819  { LocalInterpolation_ND(tk, dof2tk, Trans, I); }
1820 
1821  using FiniteElement::Project;
1822 
1823  virtual void Project(VectorCoefficient &vc,
1824  ElementTransformation &Trans, Vector &dofs) const
1825  { Project_ND(tk, dof2tk, vc, Trans, dofs); }
1826 
1827  virtual void Project(const FiniteElement &fe,
1828  ElementTransformation &Trans,
1829  DenseMatrix &I) const
1830  { Project_ND(tk, dof2tk, fe, Trans, I); }
1831 
1832  virtual void ProjectGrad(const FiniteElement &fe,
1833  ElementTransformation &Trans,
1834  DenseMatrix &grad) const
1835  { ProjectGrad_ND(tk, dof2tk, fe, Trans, grad); }
1836 
1837  virtual void ProjectCurl(const FiniteElement &fe,
1838  ElementTransformation &Trans,
1839  DenseMatrix &curl) const
1840  { ProjectCurl_ND(tk, dof2tk, fe, Trans, curl); }
1841 };
1842 
1843 
1845 {
1846  static const double tk[8];
1847 
1848  Poly_1D::Basis &cbasis1d, &obasis1d;
1849 #ifndef MFEM_THREAD_SAFE
1850  mutable Vector shape_cx, shape_ox, shape_cy, shape_oy;
1851  mutable Vector dshape_cx, dshape_cy;
1852 #endif
1853  Array<int> dof_map, dof2tk;
1854 
1855 public:
1856  ND_QuadrilateralElement(const int p);
1857  virtual void CalcVShape(const IntegrationPoint &ip,
1858  DenseMatrix &shape) const;
1859  virtual void CalcVShape(ElementTransformation &Trans,
1860  DenseMatrix &shape) const
1861  { CalcVShape_ND(Trans, shape); }
1862  virtual void CalcCurlShape(const IntegrationPoint &ip,
1863  DenseMatrix &curl_shape) const;
1865  DenseMatrix &I) const
1866  { LocalInterpolation_ND(tk, dof2tk, Trans, I); }
1867  using FiniteElement::Project;
1868  virtual void Project(VectorCoefficient &vc,
1869  ElementTransformation &Trans, Vector &dofs) const
1870  { Project_ND(tk, dof2tk, vc, Trans, dofs); }
1871  virtual void Project(const FiniteElement &fe,
1872  ElementTransformation &Trans,
1873  DenseMatrix &I) const
1874  { Project_ND(tk, dof2tk, fe, Trans, I); }
1875  virtual void ProjectGrad(const FiniteElement &fe,
1876  ElementTransformation &Trans,
1877  DenseMatrix &grad) const
1878  { ProjectGrad_ND(tk, dof2tk, fe, Trans, grad); }
1879 };
1880 
1881 
1883 {
1884  static const double tk[18], c;
1885 
1886 #ifndef MFEM_THREAD_SAFE
1887  mutable Vector shape_x, shape_y, shape_z, shape_l;
1888  mutable Vector dshape_x, dshape_y, dshape_z, dshape_l;
1889  mutable DenseMatrix u;
1890 #endif
1891  Array<int> dof2tk;
1892  DenseMatrixInverse Ti;
1893 
1894 public:
1895  ND_TetrahedronElement(const int p);
1896  virtual void CalcVShape(const IntegrationPoint &ip,
1897  DenseMatrix &shape) const;
1898  virtual void CalcVShape(ElementTransformation &Trans,
1899  DenseMatrix &shape) const
1900  { CalcVShape_ND(Trans, shape); }
1901  virtual void CalcCurlShape(const IntegrationPoint &ip,
1902  DenseMatrix &curl_shape) const;
1904  DenseMatrix &I) const
1905  { LocalInterpolation_ND(tk, dof2tk, Trans, I); }
1906  using FiniteElement::Project;
1907  virtual void Project(VectorCoefficient &vc,
1908  ElementTransformation &Trans, Vector &dofs) const
1909  { Project_ND(tk, dof2tk, vc, Trans, dofs); }
1910  virtual void Project(const FiniteElement &fe,
1911  ElementTransformation &Trans,
1912  DenseMatrix &I) const
1913  { Project_ND(tk, dof2tk, fe, Trans, I); }
1914  virtual void ProjectGrad(const FiniteElement &fe,
1915  ElementTransformation &Trans,
1916  DenseMatrix &grad) const
1917  { ProjectGrad_ND(tk, dof2tk, fe, Trans, grad); }
1918 
1919  virtual void ProjectCurl(const FiniteElement &fe,
1920  ElementTransformation &Trans,
1921  DenseMatrix &curl) const
1922  { ProjectCurl_ND(tk, dof2tk, fe, Trans, curl); }
1923 };
1924 
1926 {
1927  static const double tk[8], c;
1928 
1929 #ifndef MFEM_THREAD_SAFE
1930  mutable Vector shape_x, shape_y, shape_l;
1931  mutable Vector dshape_x, dshape_y, dshape_l;
1932  mutable DenseMatrix u;
1933  mutable Vector curlu;
1934 #endif
1935  Array<int> dof2tk;
1936  DenseMatrixInverse Ti;
1937 
1938 public:
1939  ND_TriangleElement(const int p);
1940  virtual void CalcVShape(const IntegrationPoint &ip,
1941  DenseMatrix &shape) const;
1942  virtual void CalcVShape(ElementTransformation &Trans,
1943  DenseMatrix &shape) const
1944  { CalcVShape_ND(Trans, shape); }
1945  virtual void CalcCurlShape(const IntegrationPoint &ip,
1946  DenseMatrix &curl_shape) const;
1948  DenseMatrix &I) const
1949  { LocalInterpolation_ND(tk, dof2tk, Trans, I); }
1950  using FiniteElement::Project;
1951  virtual void Project(VectorCoefficient &vc,
1952  ElementTransformation &Trans, Vector &dofs) const
1953  { Project_ND(tk, dof2tk, vc, Trans, dofs); }
1954  virtual void Project(const FiniteElement &fe,
1955  ElementTransformation &Trans,
1956  DenseMatrix &I) const
1957  { Project_ND(tk, dof2tk, fe, Trans, I); }
1958  virtual void ProjectGrad(const FiniteElement &fe,
1959  ElementTransformation &Trans,
1960  DenseMatrix &grad) const
1961  { ProjectGrad_ND(tk, dof2tk, fe, Trans, grad); }
1962 };
1963 
1964 
1966 {
1967  static const double tk[1];
1968 
1969  Poly_1D::Basis &obasis1d;
1970  Array<int> dof2tk;
1971 
1972 public:
1973  ND_SegmentElement(const int p);
1974  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
1975  { obasis1d.Eval(ip.x, shape); }
1976  virtual void CalcVShape(const IntegrationPoint &ip,
1977  DenseMatrix &shape) const;
1978  virtual void CalcVShape(ElementTransformation &Trans,
1979  DenseMatrix &shape) const
1980  { CalcVShape_ND(Trans, shape); }
1981  // virtual void CalcCurlShape(const IntegrationPoint &ip,
1982  // DenseMatrix &curl_shape) const;
1984  DenseMatrix &I) const
1985  { LocalInterpolation_ND(tk, dof2tk, Trans, I); }
1986  using FiniteElement::Project;
1987  virtual void Project(VectorCoefficient &vc,
1988  ElementTransformation &Trans, Vector &dofs) const
1989  { Project_ND(tk, dof2tk, vc, Trans, dofs); }
1990  virtual void Project(const FiniteElement &fe,
1991  ElementTransformation &Trans,
1992  DenseMatrix &I) const
1993  { Project_ND(tk, dof2tk, fe, Trans, I); }
1994  virtual void ProjectGrad(const FiniteElement &fe,
1995  ElementTransformation &Trans,
1996  DenseMatrix &grad) const
1997  { ProjectGrad_ND(tk, dof2tk, fe, Trans, grad); }
1998 };
1999 
2000 
2002 {
2003 protected:
2005  mutable int *ijk, patch, elem;
2006  mutable Vector weights;
2007 
2008 public:
2009  NURBSFiniteElement(int D, int G, int Do, int O, int F)
2010  : FiniteElement(D, G, Do, O, F)
2011  {
2012  ijk = NULL;
2013  patch = elem = -1;
2014  kv.SetSize(Dim);
2015  weights.SetSize(Dof);
2016  weights = 1.0;
2017  }
2018 
2019  void Reset () const { patch = elem = -1; }
2020  void SetIJK (int *IJK) const { ijk = IJK; }
2021  int GetPatch () const { return patch; }
2022  void SetPatch (int p) const { patch = p; }
2023  int GetElement () const { return elem; }
2024  void SetElement (int e) const { elem = e; }
2025  Array <KnotVector*> &KnotVectors() const { return kv; }
2026  Vector &Weights () const { return weights; }
2027 };
2028 
2030 {
2031 protected:
2032  mutable Vector shape_x;
2033 
2034 public:
2036  : NURBSFiniteElement(1, Geometry::SEGMENT, p + 1, p, FunctionSpace::Qk),
2037  shape_x(p + 1) { }
2038 
2039  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
2040  virtual void CalcDShape(const IntegrationPoint &ip,
2041  DenseMatrix &dshape) const;
2042 };
2043 
2045 {
2046 protected:
2048 
2049 public:
2051  : NURBSFiniteElement(2, Geometry::SQUARE, (p + 1)*(p + 1), p,
2052  FunctionSpace::Qk), u(Dof),
2053  shape_x(p + 1), shape_y(p + 1), dshape_x(p + 1), dshape_y(p + 1) { }
2054 
2055  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
2056  virtual void CalcDShape(const IntegrationPoint &ip,
2057  DenseMatrix &dshape) const;
2058 };
2059 
2061 {
2062 protected:
2064 
2065 public:
2067  : NURBSFiniteElement(3, Geometry::CUBE, (p + 1)*(p + 1)*(p + 1), p,
2068  FunctionSpace::Qk), u(Dof),
2069  shape_x(p + 1), shape_y(p + 1), shape_z(p + 1),
2070  dshape_x(p + 1), dshape_y(p + 1), dshape_z(p + 1) { }
2071 
2072  virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const;
2073  virtual void CalcDShape(const IntegrationPoint &ip,
2074  DenseMatrix &dshape) const;
2075 };
2076 
2077 }
2078 
2079 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:44
RefinedLinear3DFiniteElement()
Construct a quadratic FE on tetrahedron.
Definition: fe.cpp:4073
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1192
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1189
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7896
RT_QuadrilateralElement(const int p)
Definition: fe.cpp:8875
DenseMatrix curlshape_J
Definition: fe.hpp:280
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:5872
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1702
int GetDim() const
Returns the space dimension for the finite element.
Definition: fe.hpp:82
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:7163
L2_TetrahedronElement(const int p, const int _type=0)
Definition: fe.cpp:8564
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:2714
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1164
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3651
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:937
static void CalcBernstein(const int p, const double x, double *u, double *d)
Definition: fe.hpp:1301
Class for integration rule.
Definition: intrules.hpp:83
Linear 1D element with nodes 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:902
Quadratic 1D element with nodes the Gaussian points in [0,1] (trace of RT2)
Definition: fe.hpp:912
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:8316
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:10783
Linear3DFiniteElement()
Construct a linear FE on tetrahedron.
Definition: fe.cpp:2112
static double CalcDelta(const int p, const double x)
Definition: fe.hpp:1280
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:2896
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1732
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1776
H1_HexahedronElement(const int p)
Definition: fe.cpp:6847
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:8055
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:9646
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:30
RefinedBiLinear2DFiniteElement()
Construct a biquadratic FE on quadrilateral.
Definition: fe.cpp:4305
ND_QuadrilateralElement(const int p)
Definition: fe.cpp:10097
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1507
L2Pos_TriangleElement(const int p)
Definition: fe.cpp:8468
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.hpp:1875
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:757
H1Pos_SegmentElement(const int p)
Definition: fe.cpp:7105
int GetPatch() const
Definition: fe.hpp:2021
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.hpp:1786
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Definition: fe.cpp:116
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:649
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:4917
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:75
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:560
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:2350
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:2129
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:8552
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.hpp:1919
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:2668
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:2948
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1303
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7143
FiniteElement(D, G, Do, O, F)
L2Pos_TetrahedronElement(const int p)
Definition: fe.cpp:8715
Array< KnotVector * > kv
Definition: fe.hpp:2004
void CalcVShape_RT(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.cpp:372
void SetSize(int s)
Resizes the vector if the new size is different.
Definition: vector.hpp:259
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:9033
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:4979
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3438
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:858
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:177
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1327
void LocalInterpolation_ND(const double *tk, const Array< int > &d2t, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:687
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &h) const
Definition: fe.cpp:1685
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1596
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:2433
static void UniformPoints(const int p, double *x)
Definition: fe.cpp:6150
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
evaluate derivatives of shape function - constant 0
Definition: fe.cpp:2083
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1112
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1810
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:747
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:6675
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:10731
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7915
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1661
LagrangeHexFiniteElement(int degree)
Definition: fe.cpp:3703
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:91
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:2446
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3559
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3865
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1859
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3513
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.cpp:259
void Reset() const
Definition: fe.hpp:2019
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1741
Class for quadratic FE on triangle.
Definition: fe.hpp:486
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:10745
Basis & OpenBasis(const int p)
Definition: fe.cpp:6548
const double * ClosedPoints(const int p)
Definition: fe.cpp:6523
PositiveFiniteElement(int D, int G, int Do, int O, int F=FunctionSpace::Pk)
Definition: fe.hpp:256
Class for quadratic FE on interval.
Definition: fe.hpp:457
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1656
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:9503
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:5145
void ProjectCurl_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:509
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:2098
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:94
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1756
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:69
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1705
RefinedTriLinear3DFiniteElement()
Construct a biquadratic FE on quadrilateral.
Definition: fe.cpp:4452
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3690
Quadratic3DFiniteElement()
Construct a quadratic FE on tetrahedron.
Definition: fe.cpp:2168
int GetMapType() const
Definition: fe.hpp:98
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7769
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:5495
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:5682
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:2545
const Array< int > & GetDofMap() const
Definition: fe.hpp:1361
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:8680
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.hpp:1837
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1808
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:8770
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:2375
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:2367
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7839
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:6759
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:962
void ProjectCurl_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:547
H1Pos_HexahedronElement(const int p)
Definition: fe.cpp:7273
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.cpp:100
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:8622
RT_HexahedronElement(const int p)
Definition: fe.cpp:9083
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &h) const
Definition: fe.cpp:819
int GetElement() const
Definition: fe.hpp:2023
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7245
Linear2DFiniteElement()
Construct a linear FE on triangle.
Definition: fe.cpp:761
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1665
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:2225
Class for refined bi-linear FE on quadrilateral.
Definition: fe.hpp:1041
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:5817
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:4641
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:10628
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:728
L2_SegmentElement(const int p, const int _type=0)
Definition: fe.cpp:7740
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1828
L2Pos_HexahedronElement(const int p)
Definition: fe.cpp:8243
ND_TriangleElement(const int p)
Definition: fe.cpp:10548
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.hpp:1581
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:9682
static void CalcBasis(const int p, const double x, double *u)
Definition: fe.hpp:1264
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:784
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
evaluate shape function - constant 1
Definition: fe.cpp:2077
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &h) const
Definition: fe.cpp:1856
static void CalcBinomTerms(const int p, const double x, const double y, double *u)
Compute the terms in the expansion of the binomial (x + y)^p.
Definition: fe.cpp:6333
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &h) const
Definition: fe.cpp:63
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:3150
Cubic3DFiniteElement()
Construct a cubic FE on tetrahedron.
Definition: fe.cpp:1903
Linear1DFiniteElement()
Construct a linear FE on interval.
Definition: fe.cpp:740
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:4014
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:5805
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:1376
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:108
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:994
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1698
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:2399
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1996
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3456
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:452
Class for bilinear FE on quad with nodes at the 4 Gaussian points.
Definition: fe.hpp:432
Class for refined linear FE on interval.
Definition: fe.hpp:988
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7775
Class for refined linear FE on triangle.
Definition: fe.hpp:1008
Class for refined linear FE on tetrahedron.
Definition: fe.hpp:1028
void SetPatch(int p) const
Definition: fe.hpp:2022
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:954
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:965
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1213
Class for constant FE on triangle.
Definition: fe.hpp:625
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.hpp:1994
ND_SegmentElement(const int p)
Definition: fe.cpp:10707
H1_TriangleElement(const int p)
Definition: fe.cpp:7451
RefinedLinear1DFiniteElement()
Construct a quadratic FE on interval.
Definition: fe.cpp:3903
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:7020
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:7266
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1480
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:8410
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7518
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.hpp:1708
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3930
VectorFiniteElement(int D, int G, int Do, int O, int M, int F=FunctionSpace::Pk)
Definition: fe.hpp:331
void LocalInterpolation_RT(const double *nk, const Array< int > &d2n, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:651
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:3089
static void ChebyshevPoints(const int p, double *x)
Definition: fe.cpp:6301
NURBSFiniteElement(int D, int G, int Do, int O, int F)
Definition: fe.hpp:2009
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:5942
const Array< int > & GetDofMap() const
Definition: fe.hpp:1342
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3675
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:5783
Class for refined trilinear FE on a hexahedron.
Definition: fe.hpp:1061
ND_TetrahedronElement(const int p)
Definition: fe.cpp:10286
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:8516
Class for quadratic FE on tetrahedron.
Definition: fe.hpp:680
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1898
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:734
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1771
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1771
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1903
const IntegrationRule & GetNodes() const
Definition: fe.hpp:113
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:5032
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1035
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:8151
DenseMatrix vshape
Definition: fe.hpp:280
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:2203
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.hpp:1671
L2Pos_QuadrilateralElement(const int p)
Definition: fe.cpp:7991
P0SegmentFiniteElement(int Ord=0)
Definition: fe.cpp:2338
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:2691
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:10818
DenseMatrix curlshape
Definition: fe.hpp:280
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1907
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:864
Class for linear FE on tetrahedron.
Definition: fe.hpp:655
Class for linear FE on interval.
Definition: fe.hpp:354
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:9244
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:5083
Crouzeix-Raviart finite element on quadrilateral.
Definition: fe.hpp:728
Class for linear FE on triangle.
Definition: fe.hpp:374
Class for quadratic FE on triangle with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:510
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:8388
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:8986
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1947
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1635
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:3348
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1864
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:810
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:2138
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:9318
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7681
virtual void CalcHessian(const IntegrationPoint &ip, DenseMatrix &h) const
Definition: fe.cpp:1059
void SetElement(int e) const
Definition: fe.hpp:2024
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1137
void ProjectGrad_ND(const double *tk, const Array< int > &d2t, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.cpp:630
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:8490
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1783
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:811
virtual void ProjectCurl(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.hpp:1508
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:10842
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:1151
int GetGeomType() const
Returns the geometry type:
Definition: fe.hpp:85
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:10232
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:8016
ND_HexahedronElement(const int p)
Definition: fe.cpp:9728
IntegrationRule Nodes
Definition: fe.hpp:48
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1990
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1021
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:8271
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7123
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:7444
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:3400
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:10662
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:2560
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:4328
virtual ~FiniteElement()
Definition: fe.hpp:196
H1_QuadrilateralElement(const int p)
Definition: fe.cpp:6703
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:2457
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Overrides the scalar CalcShape function to print an error.
Definition: fe.hpp:1974
void CalcVShape_ND(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.cpp:385
DenseMatrix Jinv
Definition: fe.hpp:279
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1136
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1987
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:7786
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7706
Class for bi-quadratic FE on quadrilateral.
Definition: fe.hpp:526
Class for tri-linear FE on cube.
Definition: fe.hpp:693
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:905
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:951
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:5113
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:8738
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:723
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:772
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:4498
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:4195
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:637
static const int * Binom(const int p)
Definition: fe.cpp:6133
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:3250
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:44
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:840
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1668
static void GaussLobattoPoints(const int p, double *x)
Definition: fe.cpp:6225
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1090
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:8647
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:9944
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:2624
FiniteElement(int D, int G, int Do, int O, int F=FunctionSpace::Pk)
Definition: fe.cpp:22
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:1262
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:2283
Class for bilinear FE on quadrilateral.
Definition: fe.hpp:396
Quad1DFiniteElement()
Construct a quadratic FE on interval.
Definition: fe.cpp:946
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1978
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7399
virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const
Definition: fe.cpp:58
Class for linear FE on triangle with nodes at the 3 &quot;Gaussian&quot; points.
Definition: fe.hpp:421
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:894
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1817
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:885
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:801
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:51
H1Pos_QuadrilateralElement(const int p)
Definition: fe.cpp:7170
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1951
Array< KnotVector * > & KnotVectors() const
Definition: fe.hpp:2025
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:754
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:8439
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7420
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7845
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:230
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:2104
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
Definition: fe.hpp:565
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:838
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:8292
void Project_RT(const double *nk, const Array< int > &d2n, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:404
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:6975
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:929
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:3064
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:10434
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:2824
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:10183
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:391
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.hpp:1958
virtual void ProjectDiv(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &div) const
Definition: fe.cpp:288
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:2572
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:7856
RefinedLinear2DFiniteElement()
Construct a quadratic FE on triangle.
Definition: fe.cpp:3949
void NodalLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I, const NodalFiniteElement &fine_fe) const
Definition: fe.cpp:125
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:5737
void ProjectCurl_2D(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &curl) const
Definition: fe.cpp:158
void SetIJK(int *IJK) const
Definition: fe.hpp:2020
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1868
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.hpp:1832
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3911
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3639
RT_TetrahedronElement(const int p)
Definition: fe.cpp:9545
L2_TriangleElement(const int p, const int _type=0)
Definition: fe.cpp:8335
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.cpp:5293
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3844
Vector & Weights() const
Definition: fe.hpp:2026
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:6656
H1_TetrahedronElement(const int p)
Definition: fe.cpp:7570
BiQuad2DFiniteElement()
Construct a biquadratic FE on quadrilateral.
Definition: fe.cpp:1166
virtual void GetFaceDofs(int face, int **dofs, int *ndofs) const
Definition: fe.cpp:2158
Tensor products of 1D FEs (only degree 2 is functional)
Definition: fe.hpp:967
Basis & ClosedBasis(const int p)
Definition: fe.cpp:6571
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:2861
P0TriangleFiniteElement()
Construct P0 triangle finite element.
Definition: fe.cpp:2070
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1954
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:2344
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1744
RT_TriangleElement(const int p)
Definition: fe.cpp:9393
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:6637
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:10722
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1823
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:10018
void Eval(const double x, Vector &u) const
Definition: fe.cpp:6023
L2Pos_SegmentElement(const int p)
Definition: fe.cpp:7818
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:4108
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:5348
Class for cubic FE on tetrahedron.
Definition: fe.hpp:612
Vector data type.
Definition: vector.hpp:33
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:5279
static void CalcBernstein(const int p, const double x, double *u)
Definition: fe.hpp:1299
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1827
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:7225
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:8127
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:1087
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:9470
Quad2DFiniteElement()
Construct a quadratic FE on triangle.
Definition: fe.cpp:1004
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1871
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.hpp:1780
Describes the space on each element.
Definition: fe.hpp:27
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:2299
const double * OpenPoints(const int p)
Definition: fe.cpp:6498
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1983
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:5198
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:8035
static void GaussPoints(const int p, double *x)
Definition: fe.cpp:6165
virtual void CalcCurlShape(const IntegrationPoint &ip, DenseMatrix &curl_shape) const
Definition: fe.cpp:10476
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3696
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:1405
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:780
H1_SegmentElement(const int p)
Definition: fe.cpp:6618
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:1968
NodalFiniteElement(int D, int G, int Do, int O, int F=FunctionSpace::Pk)
Definition: fe.hpp:215
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:2410
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:6798
L2_HexahedronElement(const int p, const int _type=0)
Definition: fe.cpp:8070
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:711
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:5249
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:3468
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:415
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:4379
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:10764
virtual void GetLocalInterpolation(ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1737
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.hpp:1914
void Project_ND(const double *tk, const Array< int > &d2t, VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:566
int GetRangeType() const
Definition: fe.hpp:96
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:883
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:2766
static void CalcBasis(const int p, const double x, double *u, double *d)
Definition: fe.hpp:1273
virtual void Project(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Definition: fe.hpp:1910
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3669
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:6996
virtual void ProjectGrad(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.hpp:1747
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:6778
Basis(const int p, const double *nodes, const int _mode=1)
Definition: fe.cpp:5977
static void CalcDBinomTerms(const int p, const double x, const double y, double *d)
Definition: fe.cpp:6397
virtual void CalcVShape(const IntegrationPoint &ip, DenseMatrix &shape) const
Definition: fe.cpp:4862
TriLinear3DFiniteElement()
Construct a tri-linear FE on cube.
Definition: fe.cpp:2247
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:87
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:984
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:8106
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1693
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3966
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:3429
BiLinear2DFiniteElement()
Construct a bilinear FE on quadrilateral.
Definition: fe.cpp:788
virtual void CalcVShape(ElementTransformation &Trans, DenseMatrix &shape) const
Definition: fe.hpp:1942
Poly_1D poly1d
Definition: fe.cpp:6614
Crouzeix-Raviart finite element on triangle.
Definition: fe.hpp:716
L2_QuadrilateralElement(const int p, const int _type=0)
Definition: fe.cpp:7863
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const
Definition: fe.cpp:5925
Lagrange1DFiniteElement(int degree)
Definition: fe.cpp:3481
virtual void Project(VectorCoefficient &vc, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:2509
virtual void CalcDivShape(const IntegrationPoint &ip, Vector &divshape) const
Definition: fe.cpp:5618
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:8859
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.cpp:7935
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:850
virtual void ProjectDelta(int vertex, Vector &dofs) const
Definition: fe.hpp:673
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const
Definition: fe.cpp:7540
void ProjectGrad_RT(const double *nk, const Array< int > &d2n, const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &grad) const
Definition: fe.cpp:482
virtual void Project(Coefficient &coeff, ElementTransformation &Trans, Vector &dofs) const
Definition: fe.cpp:320