MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_coll.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_COLLECTION
13 #define MFEM_FE_COLLECTION
14 
15 #include "../config/config.hpp"
16 #include "geom.hpp"
17 #include "fe.hpp"
18 
19 namespace mfem
20 {
21 
27 {
28 public:
29  virtual const FiniteElement *
30  FiniteElementForGeometry(int GeomType) const = 0;
31 
32  virtual int DofForGeometry(int GeomType) const = 0;
33 
34  virtual int * DofOrderForOrientation(int GeomType, int Or) const = 0;
35 
36  virtual const char * Name() const { return "Undefined"; }
37 
38  int HasFaceDofs(int GeomType) const;
39 
41  int GeomType) const
42  {
43  return FiniteElementForGeometry(GeomType);
44  }
45 
47 
49 
50  static FiniteElementCollection *New(const char *name);
51 };
52 
55 {
56 protected:
57  char h1_name[32];
60  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
61 
62 public:
63  explicit H1_FECollection(const int p, const int dim = 3, const int type = 0);
64 
65  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
66  { return H1_Elements[GeomType]; }
67  virtual int DofForGeometry(int GeomType) const
68  { return H1_dof[GeomType]; }
69  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
70  virtual const char *Name() const { return h1_name; }
72 
73  virtual ~H1_FECollection();
74 };
75 
79 {
80 public:
81  explicit H1Pos_FECollection(const int p, const int dim = 3)
82  : H1_FECollection(p, dim, 1) { }
83 };
84 
89 {
90 public:
91  H1_Trace_FECollection(const int p, const int dim, const int type = 0);
92 };
93 
96 {
97 private:
98  char d_name[32];
99  FiniteElement *L2_Elements[Geometry::NumGeom];
100  FiniteElement *Tr_Elements[Geometry::NumGeom];
101  int *SegDofOrd[2]; // for rotating segment dofs in 1D
102  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
103 
104 public:
105  L2_FECollection(const int p, const int dim, const int type = 0);
106 
107  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
108  { return L2_Elements[GeomType]; }
109  virtual int DofForGeometry(int GeomType) const
110  {
111  if (L2_Elements[GeomType])
112  {
113  return L2_Elements[GeomType]->GetDof();
114  }
115  return 0;
116  }
117  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
118  virtual const char *Name() const { return d_name; }
119 
121  int GeomType) const
122  {
123  return Tr_Elements[GeomType];
124  }
125 
126  virtual ~L2_FECollection();
127 };
128 
129 // Declare an alternative name for L2_FECollection = DG_FECollection
131 
134 {
135 protected:
136  char rt_name[32];
139  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
140 
141  // Initialize only the face elements
142  void InitFaces(const int p, const int dim, const int map_type,
143  const bool signs);
144 
145  // Constructor used by the constructor of RT_Trace_FECollection
146  RT_FECollection(const int p, const int dim, const int map_type,
147  const bool signs)
148  { InitFaces(p, dim, map_type, signs); }
149 
150 public:
151  RT_FECollection(const int p, const int dim);
152 
153  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
154  { return RT_Elements[GeomType]; }
155  virtual int DofForGeometry(int GeomType) const
156  { return RT_dof[GeomType]; }
157  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
158  virtual const char *Name() const { return rt_name; }
160 
161  virtual ~RT_FECollection();
162 };
163 
168 {
169 public:
170  RT_Trace_FECollection(const int p, const int dim,
171  const int map_type = FiniteElement::INTEGRAL);
172 };
173 
178 {
179 public:
180  DG_Interface_FECollection(const int p, const int dim,
181  const int map_type = FiniteElement::VALUE);
182 };
183 
186 {
187 protected:
188  char nd_name[32];
191  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
192 
193 public:
194  ND_FECollection(const int p, const int dim);
195 
196  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
197  { return ND_Elements[GeomType]; }
198  virtual int DofForGeometry(int GeomType) const
199  { return ND_dof[GeomType]; }
200  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
201  virtual const char *Name() const { return nd_name; }
203 
204  virtual ~ND_FECollection();
205 };
206 
211 {
212 public:
213  ND_Trace_FECollection(const int p, const int dim);
214 };
215 
218 {
219 private:
220  NURBS1DFiniteElement *SegmentFE;
221  NURBS2DFiniteElement *QuadrilateralFE;
222  NURBS3DFiniteElement *ParallelepipedFE;
223 
224  char name[16];
225 
226  void Allocate(int Order);
227  void Deallocate();
228 
229 public:
230  explicit NURBSFECollection(int Order) { Allocate(Order); }
231 
232  int GetOrder() const { return SegmentFE->GetOrder(); }
233 
235  void UpdateOrder(int Order) { Deallocate(); Allocate(Order); }
236 
237  void Reset() const
238  {
239  SegmentFE->Reset();
240  QuadrilateralFE->Reset();
241  ParallelepipedFE->Reset();
242  }
243 
244  virtual const FiniteElement *
245  FiniteElementForGeometry(int GeomType) const;
246 
247  virtual int DofForGeometry(int GeomType) const;
248 
249  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
250 
251  virtual const char *Name() const { return name; }
252 
254 
255  virtual ~NURBSFECollection() { Deallocate(); }
256 };
257 
258 
261 {
262 private:
263  const PointFiniteElement PointFE;
264  const Linear1DFiniteElement SegmentFE;
265  const Linear2DFiniteElement TriangleFE;
266  const BiLinear2DFiniteElement QuadrilateralFE;
267  const Linear3DFiniteElement TetrahedronFE;
268  const TriLinear3DFiniteElement ParallelepipedFE;
269 public:
271 
272  virtual const FiniteElement *
273  FiniteElementForGeometry(int GeomType) const;
274 
275  virtual int DofForGeometry(int GeomType) const;
276 
277  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
278 
279  virtual const char * Name() const { return "Linear"; }
280 };
281 
284 {
285 private:
286  const PointFiniteElement PointFE;
287  const Quad1DFiniteElement SegmentFE;
288  const Quad2DFiniteElement TriangleFE;
289  const BiQuad2DFiniteElement QuadrilateralFE;
290  const Quadratic3DFiniteElement TetrahedronFE;
291  const LagrangeHexFiniteElement ParallelepipedFE;
292 
293 public:
294  QuadraticFECollection() : ParallelepipedFE(2) { }
295 
296  virtual const FiniteElement *
297  FiniteElementForGeometry(int GeomType) const;
298 
299  virtual int DofForGeometry(int GeomType) const;
300 
301  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
302 
303  virtual const char * Name() const { return "Quadratic"; }
304 };
305 
308 {
309 private:
310  const QuadPos1DFiniteElement SegmentFE;
311  const BiQuadPos2DFiniteElement QuadrilateralFE;
312 
313 public:
315 
316  virtual const FiniteElement *
317  FiniteElementForGeometry(int GeomType) const;
318 
319  virtual int DofForGeometry(int GeomType) const;
320 
321  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
322 
323  virtual const char * Name() const { return "QuadraticPos"; }
324 };
325 
328 {
329 private:
330  const PointFiniteElement PointFE;
331  const Cubic1DFiniteElement SegmentFE;
332  const Cubic2DFiniteElement TriangleFE;
333  const BiCubic2DFiniteElement QuadrilateralFE;
334  const Cubic3DFiniteElement TetrahedronFE;
335  const LagrangeHexFiniteElement ParallelepipedFE;
336 
337 public:
338  CubicFECollection() : ParallelepipedFE(3) { }
339 
340  virtual const FiniteElement *
341  FiniteElementForGeometry(int GeomType) const;
342 
343  virtual int DofForGeometry(int GeomType) const;
344 
345  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
346 
347  virtual const char * Name() const { return "Cubic"; }
348 };
349 
352 {
353 private:
354  const P0SegmentFiniteElement SegmentFE;
355  const CrouzeixRaviartFiniteElement TriangleFE;
356  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
357 public:
358  CrouzeixRaviartFECollection() : SegmentFE(1) { }
359 
360  virtual const FiniteElement *
361  FiniteElementForGeometry(int GeomType) const;
362 
363  virtual int DofForGeometry(int GeomType) const;
364 
365  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
366 
367  virtual const char * Name() const { return "CrouzeixRaviart"; }
368 };
369 
372 {
373 private:
374  const P0TriangleFiniteElement TriangleFE;
375  const P1TetNonConfFiniteElement TetrahedronFE;
376  const P0QuadFiniteElement QuadrilateralFE;
377  const RotTriLinearHexFiniteElement ParallelepipedFE;
378 
379 public:
381 
382  virtual const FiniteElement *
383  FiniteElementForGeometry(int GeomType) const;
384 
385  virtual int DofForGeometry(int GeomType) const;
386 
387  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
388 
389  virtual const char * Name() const { return "LinearNonConf3D"; }
390 };
391 
392 
396 {
397 private:
398  const P0SegmentFiniteElement SegmentFE; // normal component on edge
399  const RT0TriangleFiniteElement TriangleFE;
400  const RT0QuadFiniteElement QuadrilateralFE;
401 public:
402  RT0_2DFECollection() : SegmentFE(0) { }
403 
404  virtual const FiniteElement *
405  FiniteElementForGeometry(int GeomType) const;
406 
407  virtual int DofForGeometry(int GeomType) const;
408 
409  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
410 
411  virtual const char * Name() const { return "RT0_2D"; }
412 };
413 
417 {
418 private:
419  const P1SegmentFiniteElement SegmentFE; // normal component on edge
420  const RT1TriangleFiniteElement TriangleFE;
421  const RT1QuadFiniteElement QuadrilateralFE;
422 public:
424 
425  virtual const FiniteElement *
426  FiniteElementForGeometry(int GeomType) const;
427 
428  virtual int DofForGeometry(int GeomType) const;
429 
430  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
431 
432  virtual const char * Name() const { return "RT1_2D"; }
433 };
434 
438 {
439 private:
440  const P2SegmentFiniteElement SegmentFE; // normal component on edge
441  const RT2TriangleFiniteElement TriangleFE;
442  const RT2QuadFiniteElement QuadrilateralFE;
443 public:
445 
446  virtual const FiniteElement *
447  FiniteElementForGeometry(int GeomType) const;
448 
449  virtual int DofForGeometry(int GeomType) const;
450 
451  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
452 
453  virtual const char * Name() const { return "RT2_2D"; }
454 };
455 
459 {
460 private:
461  const P0TriangleFiniteElement TriangleFE;
462  const P0QuadFiniteElement QuadrilateralFE;
463 public:
465 
466  virtual const FiniteElement *
467  FiniteElementForGeometry(int GeomType) const;
468 
469  virtual int DofForGeometry(int GeomType) const;
470 
471  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
472 
473  virtual const char * Name() const { return "Const2D"; }
474 };
475 
479 {
480 private:
481  const Linear2DFiniteElement TriangleFE;
482  const BiLinear2DFiniteElement QuadrilateralFE;
483 
484 public:
486 
487  virtual const FiniteElement *
488  FiniteElementForGeometry(int GeomType) const;
489 
490  virtual int DofForGeometry(int GeomType) const;
491 
492  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
493 
494  virtual const char * Name() const { return "LinearDiscont2D"; }
495 };
496 
499 {
500 private:
501  // const CrouzeixRaviartFiniteElement TriangleFE;
502  const GaussLinear2DFiniteElement TriangleFE;
503  const GaussBiLinear2DFiniteElement QuadrilateralFE;
504 
505 public:
507 
508  virtual const FiniteElement *
509  FiniteElementForGeometry(int GeomType) const;
510 
511  virtual int DofForGeometry(int GeomType) const;
512 
513  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
514 
515  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
516 };
517 
520 {
521 private:
522  const P1OnQuadFiniteElement QuadrilateralFE;
523 public:
525  virtual const FiniteElement *
526  FiniteElementForGeometry(int GeomType) const;
527  virtual int DofForGeometry(int GeomType) const;
528  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
529  virtual const char * Name() const { return "P1OnQuad"; }
530 };
531 
535 {
536 private:
537  const Quad2DFiniteElement TriangleFE;
538  const BiQuad2DFiniteElement QuadrilateralFE;
539 
540 public:
542 
543  virtual const FiniteElement *
544  FiniteElementForGeometry(int GeomType) const;
545 
546  virtual int DofForGeometry(int GeomType) const;
547 
548  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
549 
550  virtual const char * Name() const { return "QuadraticDiscont2D"; }
551 };
552 
555 {
556 private:
557  const BiQuadPos2DFiniteElement QuadrilateralFE;
558 
559 public:
561  virtual const FiniteElement *
562  FiniteElementForGeometry(int GeomType) const;
563  virtual int DofForGeometry(int GeomType) const;
564  virtual int * DofOrderForOrientation(int GeomType, int Or) const
565  { return NULL; }
566  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
567 };
568 
571 {
572 private:
573  // const Quad2DFiniteElement TriangleFE;
574  const GaussQuad2DFiniteElement TriangleFE;
575  const GaussBiQuad2DFiniteElement QuadrilateralFE;
576 
577 public:
579 
580  virtual const FiniteElement *
581  FiniteElementForGeometry(int GeomType) const;
582 
583  virtual int DofForGeometry(int GeomType) const;
584 
585  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
586 
587  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
588 };
589 
593 {
594 private:
595  const Cubic2DFiniteElement TriangleFE;
596  const BiCubic2DFiniteElement QuadrilateralFE;
597 
598 public:
600 
601  virtual const FiniteElement *
602  FiniteElementForGeometry(int GeomType) const;
603 
604  virtual int DofForGeometry(int GeomType) const;
605 
606  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
607 
608  virtual const char * Name() const { return "CubicDiscont2D"; }
609 };
610 
614 {
615 private:
616  const P0TetFiniteElement TetrahedronFE;
617  const P0HexFiniteElement ParallelepipedFE;
618 
619 public:
621 
622  virtual const FiniteElement *
623  FiniteElementForGeometry(int GeomType) const;
624 
625  virtual int DofForGeometry(int GeomType) const;
626 
627  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
628 
629  virtual const char * Name() const { return "Const3D"; }
630 };
631 
635 {
636 private:
637  const Linear3DFiniteElement TetrahedronFE;
638  const TriLinear3DFiniteElement ParallelepipedFE;
639 
640 public:
642 
643  virtual const FiniteElement *
644  FiniteElementForGeometry(int GeomType) const;
645 
646  virtual int DofForGeometry(int GeomType) const;
647 
648  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
649 
650  virtual const char * Name() const { return "LinearDiscont3D"; }
651 };
652 
656 {
657 private:
658  const Quadratic3DFiniteElement TetrahedronFE;
659  const LagrangeHexFiniteElement ParallelepipedFE;
660 
661 public:
662  QuadraticDiscont3DFECollection () : ParallelepipedFE(2) { }
663 
664  virtual const FiniteElement *
665  FiniteElementForGeometry(int GeomType) const;
666 
667  virtual int DofForGeometry(int GeomType) const;
668 
669  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
670 
671  virtual const char * Name() const { return "QuadraticDiscont3D"; }
672 };
673 
676 {
677 private:
678  const PointFiniteElement PointFE;
679  const RefinedLinear1DFiniteElement SegmentFE;
680  const RefinedLinear2DFiniteElement TriangleFE;
681  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
682  const RefinedLinear3DFiniteElement TetrahedronFE;
683  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
684 
685 public:
687 
688  virtual const FiniteElement *
689  FiniteElementForGeometry(int GeomType) const;
690 
691  virtual int DofForGeometry(int GeomType) const;
692 
693  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
694 
695  virtual const char * Name() const { return "RefinedLinear"; }
696 };
697 
701 {
702 private:
703  const Nedelec1HexFiniteElement HexahedronFE;
704  const Nedelec1TetFiniteElement TetrahedronFE;
705 
706 public:
708 
709  virtual const FiniteElement *
710  FiniteElementForGeometry(int GeomType) const;
711 
712  virtual int DofForGeometry(int GeomType) const;
713 
714  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
715 
716  virtual const char * Name() const { return "ND1_3D"; }
717 };
718 
722 {
723 private:
724  const P0TriangleFiniteElement TriangleFE;
725  const P0QuadFiniteElement QuadrilateralFE;
726  const RT0HexFiniteElement HexahedronFE;
727  const RT0TetFiniteElement TetrahedronFE;
728 public:
730 
731  virtual const FiniteElement *
732  FiniteElementForGeometry(int GeomType) const;
733 
734  virtual int DofForGeometry(int GeomType) const;
735 
736  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
737 
738  virtual const char * Name() const { return "RT0_3D"; }
739 };
740 
744 {
745 private:
746  const Linear2DFiniteElement TriangleFE;
747  const BiLinear2DFiniteElement QuadrilateralFE;
748  const RT1HexFiniteElement HexahedronFE;
749 public:
751 
752  virtual const FiniteElement *
753  FiniteElementForGeometry(int GeomType) const;
754 
755  virtual int DofForGeometry(int GeomType) const;
756 
757  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
758 
759  virtual const char * Name() const { return "RT1_3D"; }
760 };
761 
764 {
765 private:
766  char d_name[32];
767  int GeomType;
768  FiniteElement *Local_Element;
769 
770 public:
771  Local_FECollection(const char *fe_name);
772 
773  virtual const FiniteElement *FiniteElementForGeometry(int _GeomType) const
774  { return (GeomType == _GeomType) ? Local_Element : NULL; }
775  virtual int DofForGeometry(int _GeomType) const
776  { return (GeomType == _GeomType) ? Local_Element->GetDof() : 0; }
777  virtual int *DofOrderForOrientation(int GeomType, int Or) const
778  { return NULL; }
779  virtual const char *Name() const { return d_name; }
780 
781  virtual ~Local_FECollection() { delete Local_Element; }
782 };
783 
784 }
785 
786 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:44
int GetOrder() const
Definition: fe_coll.hpp:232
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:217
virtual const char * Name() const
Definition: fe_coll.hpp:347
virtual const char * Name() const
Definition: fe_coll.hpp:671
virtual const char * Name() const
Definition: fe_coll.hpp:587
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:929
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:570
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:138
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:189
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1034
void UpdateOrder(int Order)
Change the order of the collection.
Definition: fe_coll.hpp:235
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 const char * Name() const
Definition: fe_coll.hpp:695
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:38
virtual const char * Name() const
Definition: fe_coll.hpp:629
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:301
virtual const char * Name() const
Definition: fe_coll.hpp:494
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1938
static const int NumGeom
Definition: geom.hpp:34
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1932
virtual const char * Name() const
Definition: fe_coll.hpp:303
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1097
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:857
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:775
virtual ~L2_FECollection()
Definition: fe_coll.cpp:1460
virtual const char * Name() const
Definition: fe_coll.hpp:529
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:390
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1018
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:620
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:678
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:656
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:586
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:91
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:260
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1125
void Reset() const
Definition: fe.hpp:2019
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:957
ND_FECollection(const int p, const int dim)
Definition: fe_coll.cpp:1679
Class for quadratic FE on triangle.
Definition: fe.hpp:486
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:496
Class for quadratic FE on interval.
Definition: fe.hpp:457
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:964
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1063
virtual const char * Name() const
Definition: fe_coll.hpp:158
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1834
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:155
L2_FECollection(const int p, const int dim, const int type=0)
Definition: fe_coll.cpp:1338
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:809
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:795
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:578
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:1502
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1631
virtual const char * Name() const
Definition: fe_coll.hpp:251
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:996
virtual const char * Name() const
Definition: fe_coll.hpp:411
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:893
virtual const char * Name() const
Definition: fe_coll.hpp:608
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:523
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:746
Finite element collection on a macro-element.
Definition: fe_coll.hpp:675
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:469
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:225
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:248
Class for refined bi-linear FE on quadrilateral.
Definition: fe.hpp:1041
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:196
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:510
virtual const char * Name() const
Definition: fe_coll.hpp:201
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:648
virtual ~RT_FECollection()
Definition: fe_coll.cpp:1636
virtual int DofForGeometry(int _GeomType) const
Definition: fe_coll.hpp:775
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:564
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1610
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1079
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:327
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:612
virtual const char * Name() const
Definition: fe_coll.hpp:432
int dim
Definition: ex3.cpp:48
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 int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1110
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:59
Class for refined linear FE on triangle.
Definition: fe.hpp:1008
Class for refined linear FE on tetrahedron.
Definition: fe.hpp:1028
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:137
Class for constant FE on triangle.
Definition: fe.hpp:625
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:314
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:598
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:241
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:665
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:869
ND_Trace_FECollection(const int p, const int dim)
Definition: fe_coll.cpp:1852
Class for refined trilinear FE on a hexahedron.
Definition: fe.hpp:1061
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:763
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:289
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:554
virtual ~ND_FECollection()
Definition: fe_coll.cpp:1840
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:885
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:427
Class for quadratic FE on tetrahedron.
Definition: fe.hpp:680
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:777
RT_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL)
Definition: fe_coll.cpp:1647
H1_FECollection(const int p, const int dim=3, const int type=0)
Definition: fe_coll.cpp:1147
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:698
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:552
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:307
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:817
virtual ~Local_FECollection()
Definition: fe_coll.hpp:781
virtual const char * Name() const
Definition: fe_coll.hpp:515
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:980
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1311
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:209
Class for linear FE on tetrahedron.
Definition: fe.hpp:655
Class for linear FE on interval.
Definition: fe.hpp:354
Crouzeix-Raviart finite element on quadrilateral.
Definition: fe.hpp:728
Class for linear FE on triangle.
Definition: fe.hpp:374
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1049
Class for quadratic FE on triangle with nodes at the "Gaussian" points.
Definition: fe.hpp:510
virtual const char * Name() const
Definition: fe_coll.hpp:70
virtual const char * Name() const
Definition: fe_coll.hpp:650
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1297
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:58
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:371
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:356
DG_Interface_FECollection(const int p, const int dim, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1663
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:351
NURBSFECollection(int Order)
Definition: fe_coll.hpp:230
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:133
L2_FECollection DG_FECollection
Definition: fe_coll.hpp:130
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:417
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:81
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:403
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:482
virtual const char * Name() const
Definition: fe_coll.hpp:367
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:109
virtual const char * Name() const
Definition: fe_coll.hpp:279
virtual const char * Name() const
Definition: fe_coll.hpp:453
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:120
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:905
virtual const char * Name() const
Definition: fe_coll.hpp:389
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:712
virtual const char * Name() const
Definition: fe_coll.hpp:716
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:88
Class for bi-quadratic FE on quadrilateral.
Definition: fe.hpp:526
Class for tri-linear FE on cube.
Definition: fe.hpp:693
H1_Trace_FECollection(const int p, const int dim, const int type=0)
Definition: fe_coll.cpp:1323
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:720
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:280
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:633
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:153
virtual const char * Name() const
Definition: fe_coll.hpp:550
Class for bilinear FE on quadrilateral.
Definition: fe.hpp:396
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:941
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:783
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:760
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:498
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1913
virtual int * DofOrderForOrientation(int GeomType, int Or) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:118
static FiniteElementCollection * New(const char *name)
Definition: fe_coll.cpp:44
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:40
Class for linear FE on triangle with nodes at the 3 "Gaussian" points.
Definition: fe.hpp:421
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:921
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1006
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:198
virtual const char * Name() const
Definition: fe_coll.hpp:36
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:264
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:65
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
Definition: fe.hpp:565
virtual const char * Name() const
Definition: fe_coll.hpp:473
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:537
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1807
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1276
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:454
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:67
RT_FECollection(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.hpp:146
virtual const char * Name() const
Definition: fe_coll.hpp:738
Tensor products of 1D FEs (only degree 2 is functional)
Definition: fe.hpp:967
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:283
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:107
int HasFaceDofs(int GeomType) const
Definition: fe_coll.cpp:25
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:185
void Reset() const
Definition: fe_coll.hpp:237
Class for cubic FE on tetrahedron.
Definition: fe.hpp:612
virtual const char * Name() const
Definition: fe_coll.hpp:779
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:847
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1443
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:54
virtual const char * Name() const
Definition: fe_coll.hpp:759
virtual const char * Name() const
Definition: fe_coll.hpp:566
virtual const FiniteElement * FiniteElementForGeometry(int _GeomType) const
Definition: fe_coll.hpp:773
virtual int DofForGeometry(int GeomType) const =0
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:564
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1926
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:1859
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:731
virtual const char * Name() const
Definition: fe_coll.hpp:323
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:831
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:519
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:340
virtual ~FiniteElementCollection()
Definition: fe_coll.hpp:48
virtual ~NURBSFECollection()
Definition: fe_coll.hpp:255
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:190
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:440
Crouzeix-Raviart finite element on triangle.
Definition: fe.hpp:716
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:324
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const =0
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:95
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:686