MFEM  v3.0
 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.googlecode.com.
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 
48  static FiniteElementCollection *New(const char *name);
49 };
50 
53 {
54 private:
55  char h1_name[32];
56  FiniteElement *H1_Elements[Geometry::NumGeom];
57  int H1_dof[Geometry::NumGeom];
58  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
59 
60 public:
61  explicit H1_FECollection(const int p, const int dim = 3, const int type = 0);
62 
63  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
64  { return H1_Elements[GeomType]; }
65  virtual int DofForGeometry(int GeomType) const
66  { return H1_dof[GeomType]; }
67  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
68  virtual const char *Name() const { return h1_name; }
69 
70  virtual ~H1_FECollection();
71 };
72 
76 {
77 public:
78  explicit H1Pos_FECollection(const int p, const int dim = 3)
79  : H1_FECollection(p, dim, 1) { }
80 };
81 
84 {
85 private:
86  char d_name[32];
87  FiniteElement *L2_Elements[Geometry::NumGeom];
88  FiniteElement *Tr_Elements[Geometry::NumGeom];
89  int *SegDofOrd[2]; // for rotating segment dofs in 1D
90  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
91 
92 public:
93  L2_FECollection(const int p, const int dim, const int type = 0);
94 
95  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
96  { return L2_Elements[GeomType]; }
97  virtual int DofForGeometry(int GeomType) const
98  {
99  if (L2_Elements[GeomType])
100  return L2_Elements[GeomType]->GetDof();
101  return 0;
102  }
103  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
104  virtual const char *Name() const { return d_name; }
105 
107  int GeomType) const
108  {
109  return Tr_Elements[GeomType];
110  }
111 
112  virtual ~L2_FECollection();
113 };
114 
115 // Declare an alternative name for L2_FECollection = DG_FECollection
117 
120 {
121 protected:
122  char rt_name[32];
125  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
126 
127  // Initialize only the face elements
128  void InitFaces(const int p, const int dim, const int map_type);
129 
130  // Constructor used by the constructor of RT_Trace_FECollection
131  RT_FECollection(const int p, const int dim, const int map_type)
132  { InitFaces(p, dim, map_type); }
133 
134 public:
135  RT_FECollection(const int p, const int dim);
136 
137  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
138  { return RT_Elements[GeomType]; }
139  virtual int DofForGeometry(int GeomType) const
140  { return RT_dof[GeomType]; }
141  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
142  virtual const char *Name() const { return rt_name; }
143 
144  virtual ~RT_FECollection();
145 };
146 
151 {
152 public:
153  RT_Trace_FECollection(const int p, const int dim,
154  const int map_type = FiniteElement::INTEGRAL);
155 };
156 
159 {
160 private:
161  char nd_name[32];
162  FiniteElement *ND_Elements[Geometry::NumGeom];
163  int ND_dof[Geometry::NumGeom];
164  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
165 
166 public:
167  ND_FECollection(const int p, const int dim);
168 
169  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
170  { return ND_Elements[GeomType]; }
171  virtual int DofForGeometry(int GeomType) const
172  { return ND_dof[GeomType]; }
173  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
174  virtual const char *Name() const { return nd_name; }
175 
176  virtual ~ND_FECollection();
177 };
178 
181 {
182 private:
183  NURBS1DFiniteElement *SegmentFE;
184  NURBS2DFiniteElement *QuadrilateralFE;
185  NURBS3DFiniteElement *ParallelepipedFE;
186 
187  char name[16];
188 
189  void Allocate(int Order);
190  void Deallocate();
191 
192 public:
193  explicit NURBSFECollection(int Order) { Allocate(Order); }
194 
195  int GetOrder() const { return SegmentFE->GetOrder(); }
196 
198  void UpdateOrder(int Order) { Deallocate(); Allocate(Order); }
199 
200  void Reset() const
201  {
202  SegmentFE->Reset();
203  QuadrilateralFE->Reset();
204  ParallelepipedFE->Reset();
205  }
206 
207  virtual const FiniteElement *
208  FiniteElementForGeometry(int GeomType) const;
209 
210  virtual int DofForGeometry(int GeomType) const;
211 
212  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
213 
214  virtual const char *Name() const { return name; }
215 
216  virtual ~NURBSFECollection() { Deallocate(); }
217 };
218 
219 
222 {
223 private:
224  const PointFiniteElement PointFE;
225  const Linear1DFiniteElement SegmentFE;
226  const Linear2DFiniteElement TriangleFE;
227  const BiLinear2DFiniteElement QuadrilateralFE;
228  const Linear3DFiniteElement TetrahedronFE;
229  const TriLinear3DFiniteElement ParallelepipedFE;
230 public:
232 
233  virtual const FiniteElement *
234  FiniteElementForGeometry(int GeomType) const;
235 
236  virtual int DofForGeometry(int GeomType) const;
237 
238  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
239 
240  virtual const char * Name() const { return "Linear"; };
241 };
242 
245 {
246 private:
247  const PointFiniteElement PointFE;
248  const Quad1DFiniteElement SegmentFE;
249  const Quad2DFiniteElement TriangleFE;
250  const BiQuad2DFiniteElement QuadrilateralFE;
251  const Quadratic3DFiniteElement TetrahedronFE;
252  const LagrangeHexFiniteElement ParallelepipedFE;
253 
254 public:
255  QuadraticFECollection() : ParallelepipedFE(2) { };
256 
257  virtual const FiniteElement *
258  FiniteElementForGeometry(int GeomType) const;
259 
260  virtual int DofForGeometry(int GeomType) const;
261 
262  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
263 
264  virtual const char * Name() const { return "Quadratic"; };
265 };
266 
269 {
270 private:
271  const QuadPos1DFiniteElement SegmentFE;
272  const BiQuadPos2DFiniteElement QuadrilateralFE;
273 
274 public:
276 
277  virtual const FiniteElement *
278  FiniteElementForGeometry(int GeomType) const;
279 
280  virtual int DofForGeometry(int GeomType) const;
281 
282  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
283 
284  virtual const char * Name() const { return "QuadraticPos"; };
285 };
286 
289 {
290 private:
291  const PointFiniteElement PointFE;
292  const Cubic1DFiniteElement SegmentFE;
293  const Cubic2DFiniteElement TriangleFE;
294  const BiCubic2DFiniteElement QuadrilateralFE;
295  const Cubic3DFiniteElement TetrahedronFE;
296  const LagrangeHexFiniteElement ParallelepipedFE;
297 
298 public:
299  CubicFECollection() : ParallelepipedFE(3) { };
300 
301  virtual const FiniteElement *
302  FiniteElementForGeometry(int GeomType) const;
303 
304  virtual int DofForGeometry(int GeomType) const;
305 
306  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
307 
308  virtual const char * Name() const { return "Cubic"; };
309 };
310 
313 {
314 private:
315  const P0SegmentFiniteElement SegmentFE;
316  const CrouzeixRaviartFiniteElement TriangleFE;
317  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
318 public:
319  CrouzeixRaviartFECollection() : SegmentFE(1) { };
320 
321  virtual const FiniteElement *
322  FiniteElementForGeometry(int GeomType) const;
323 
324  virtual int DofForGeometry(int GeomType) const;
325 
326  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
327 
328  virtual const char * Name() const { return "CrouzeixRaviart"; };
329 };
330 
333 {
334 private:
335  const P0TriangleFiniteElement TriangleFE;
336  const P1TetNonConfFiniteElement TetrahedronFE;
337  const P0QuadFiniteElement QuadrilateralFE;
338  const RotTriLinearHexFiniteElement ParallelepipedFE;
339 
340 public:
342 
343  virtual const FiniteElement *
344  FiniteElementForGeometry(int GeomType) const;
345 
346  virtual int DofForGeometry(int GeomType) const;
347 
348  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
349 
350  virtual const char * Name() const { return "LinearNonConf3D"; };
351 };
352 
353 
357 {
358 private:
359  const P0SegmentFiniteElement SegmentFE; // normal component on edge
360  const RT0TriangleFiniteElement TriangleFE;
361  const RT0QuadFiniteElement QuadrilateralFE;
362 public:
363  RT0_2DFECollection() : SegmentFE(0) { };
364 
365  virtual const FiniteElement *
366  FiniteElementForGeometry(int GeomType) const;
367 
368  virtual int DofForGeometry(int GeomType) const;
369 
370  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
371 
372  virtual const char * Name() const { return "RT0_2D"; };
373 };
374 
378 {
379 private:
380  const P1SegmentFiniteElement SegmentFE; // normal component on edge
381  const RT1TriangleFiniteElement TriangleFE;
382  const RT1QuadFiniteElement QuadrilateralFE;
383 public:
385 
386  virtual const FiniteElement *
387  FiniteElementForGeometry(int GeomType) const;
388 
389  virtual int DofForGeometry(int GeomType) const;
390 
391  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
392 
393  virtual const char * Name() const { return "RT1_2D"; };
394 };
395 
399 {
400 private:
401  const P2SegmentFiniteElement SegmentFE; // normal component on edge
402  const RT2TriangleFiniteElement TriangleFE;
403  const RT2QuadFiniteElement QuadrilateralFE;
404 public:
406 
407  virtual const FiniteElement *
408  FiniteElementForGeometry(int GeomType) const;
409 
410  virtual int DofForGeometry(int GeomType) const;
411 
412  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
413 
414  virtual const char * Name() const { return "RT2_2D"; };
415 };
416 
420 {
421 private:
422  const P0TriangleFiniteElement TriangleFE;
423  const P0QuadFiniteElement QuadrilateralFE;
424 public:
426 
427  virtual const FiniteElement *
428  FiniteElementForGeometry(int GeomType) const;
429 
430  virtual int DofForGeometry(int GeomType) const;
431 
432  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
433 
434  virtual const char * Name() const { return "Const2D"; };
435 };
436 
440 {
441 private:
442  const Linear2DFiniteElement TriangleFE;
443  const BiLinear2DFiniteElement QuadrilateralFE;
444 
445 public:
447 
448  virtual const FiniteElement *
449  FiniteElementForGeometry(int GeomType) const;
450 
451  virtual int DofForGeometry(int GeomType) const;
452 
453  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
454 
455  virtual const char * Name() const { return "LinearDiscont2D"; };
456 };
457 
460 {
461 private:
462  // const CrouzeixRaviartFiniteElement TriangleFE;
463  const GaussLinear2DFiniteElement TriangleFE;
464  const GaussBiLinear2DFiniteElement QuadrilateralFE;
465 
466 public:
468 
469  virtual const FiniteElement *
470  FiniteElementForGeometry(int GeomType) const;
471 
472  virtual int DofForGeometry(int GeomType) const;
473 
474  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
475 
476  virtual const char * Name() const { return "GaussLinearDiscont2D"; };
477 };
478 
481 {
482 private:
483  const P1OnQuadFiniteElement QuadrilateralFE;
484 public:
486  virtual const FiniteElement *
487  FiniteElementForGeometry(int GeomType) const;
488  virtual int DofForGeometry(int GeomType) const;
489  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
490  virtual const char * Name() const { return "P1OnQuad"; };
491 };
492 
496 {
497 private:
498  const Quad2DFiniteElement TriangleFE;
499  const BiQuad2DFiniteElement QuadrilateralFE;
500 
501 public:
503 
504  virtual const FiniteElement *
505  FiniteElementForGeometry(int GeomType) const;
506 
507  virtual int DofForGeometry(int GeomType) const;
508 
509  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
510 
511  virtual const char * Name() const { return "QuadraticDiscont2D"; };
512 };
513 
516 {
517 private:
518  const BiQuadPos2DFiniteElement QuadrilateralFE;
519 
520 public:
522  virtual const FiniteElement *
523  FiniteElementForGeometry(int GeomType) const;
524  virtual int DofForGeometry(int GeomType) const;
525  virtual int * DofOrderForOrientation(int GeomType, int Or) const
526  { return NULL; };
527  virtual const char * Name() const { return "QuadraticPosDiscont2D"; };
528 };
529 
532 {
533 private:
534  // const Quad2DFiniteElement TriangleFE;
535  const GaussQuad2DFiniteElement TriangleFE;
536  const GaussBiQuad2DFiniteElement QuadrilateralFE;
537 
538 public:
540 
541  virtual const FiniteElement *
542  FiniteElementForGeometry(int GeomType) const;
543 
544  virtual int DofForGeometry(int GeomType) const;
545 
546  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
547 
548  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; };
549 };
550 
554 {
555 private:
556  const Cubic2DFiniteElement TriangleFE;
557  const BiCubic2DFiniteElement QuadrilateralFE;
558 
559 public:
561 
562  virtual const FiniteElement *
563  FiniteElementForGeometry(int GeomType) const;
564 
565  virtual int DofForGeometry(int GeomType) const;
566 
567  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
568 
569  virtual const char * Name() const { return "CubicDiscont2D"; };
570 };
571 
575 {
576 private:
577  const P0TetFiniteElement TetrahedronFE;
578  const P0HexFiniteElement ParallelepipedFE;
579 
580 public:
582 
583  virtual const FiniteElement *
584  FiniteElementForGeometry(int GeomType) const;
585 
586  virtual int DofForGeometry(int GeomType) const;
587 
588  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
589 
590  virtual const char * Name() const { return "Const3D"; };
591 };
592 
596 {
597 private:
598  const Linear3DFiniteElement TetrahedronFE;
599  const TriLinear3DFiniteElement ParallelepipedFE;
600 
601 public:
603 
604  virtual const FiniteElement *
605  FiniteElementForGeometry(int GeomType) const;
606 
607  virtual int DofForGeometry(int GeomType) const;
608 
609  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
610 
611  virtual const char * Name() const { return "LinearDiscont3D"; };
612 };
613 
617 {
618 private:
619  const Quadratic3DFiniteElement TetrahedronFE;
620  const LagrangeHexFiniteElement ParallelepipedFE;
621 
622 public:
623  QuadraticDiscont3DFECollection () : ParallelepipedFE(2) { };
624 
625  virtual const FiniteElement *
626  FiniteElementForGeometry(int GeomType) const;
627 
628  virtual int DofForGeometry(int GeomType) const;
629 
630  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
631 
632  virtual const char * Name() const { return "QuadraticDiscont3D"; };
633 };
634 
637 {
638 private:
639  const PointFiniteElement PointFE;
640  const RefinedLinear1DFiniteElement SegmentFE;
641  const RefinedLinear2DFiniteElement TriangleFE;
642  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
643  const RefinedLinear3DFiniteElement TetrahedronFE;
644  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
645 
646 public:
648 
649  virtual const FiniteElement *
650  FiniteElementForGeometry(int GeomType) const;
651 
652  virtual int DofForGeometry(int GeomType) const;
653 
654  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
655 
656  virtual const char * Name() const { return "RefinedLinear"; };
657 };
658 
662 {
663 private:
664  const Nedelec1HexFiniteElement HexahedronFE;
665  const Nedelec1TetFiniteElement TetrahedronFE;
666 
667 public:
669 
670  virtual const FiniteElement *
671  FiniteElementForGeometry(int GeomType) const;
672 
673  virtual int DofForGeometry(int GeomType) const;
674 
675  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
676 
677  virtual const char * Name() const { return "ND1_3D"; };
678 };
679 
683 {
684 private:
685  const P0TriangleFiniteElement TriangleFE;
686  const P0QuadFiniteElement QuadrilateralFE;
687  const RT0HexFiniteElement HexahedronFE;
688  const RT0TetFiniteElement TetrahedronFE;
689 public:
691 
692  virtual const FiniteElement *
693  FiniteElementForGeometry(int GeomType) const;
694 
695  virtual int DofForGeometry(int GeomType) const;
696 
697  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
698 
699  virtual const char * Name() const { return "RT0_3D"; };
700 };
701 
705 {
706 private:
707  const Linear2DFiniteElement TriangleFE;
708  const BiLinear2DFiniteElement QuadrilateralFE;
709  const RT1HexFiniteElement HexahedronFE;
710 public:
712 
713  virtual const FiniteElement *
714  FiniteElementForGeometry(int GeomType) const;
715 
716  virtual int DofForGeometry(int GeomType) const;
717 
718  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
719 
720  virtual const char * Name() const { return "RT1_3D"; };
721 };
722 
725 {
726 private:
727  char d_name[32];
728  int GeomType;
729  FiniteElement *Local_Element;
730 
731 public:
732  Local_FECollection(const char *fe_name);
733 
734  virtual const FiniteElement *FiniteElementForGeometry(int _GeomType) const
735  { return (GeomType == _GeomType) ? Local_Element : NULL; }
736  virtual int DofForGeometry(int _GeomType) const
737  { return (GeomType == _GeomType) ? Local_Element->GetDof() : 0; }
738  virtual int *DofOrderForOrientation(int GeomType, int Or) const
739  { return NULL; }
740  virtual const char *Name() const { return d_name; }
741 
742  virtual ~Local_FECollection() { delete Local_Element; }
743 };
744 
745 }
746 
747 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:42
int GetOrder() const
Definition: fe_coll.hpp:195
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:180
virtual const char * Name() const
Definition: fe_coll.hpp:308
virtual const char * Name() const
Definition: fe_coll.hpp:632
virtual const char * Name() const
Definition: fe_coll.hpp:548
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:829
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:531
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:124
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:933
void UpdateOrder(int Order)
Change the order of the collection.
Definition: fe_coll.hpp:198
Linear 1D element with nodes 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:890
Quadratic 1D element with nodes the Gaussian points in [0,1] (trace of RT2)
Definition: fe.hpp:900
virtual const char * Name() const
Definition: fe_coll.hpp:656
virtual const char * Name() const
Definition: fe_coll.hpp:590
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:210
virtual const char * Name() const
Definition: fe_coll.hpp:455
static const int NumGeom
Definition: geom.hpp:34
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1673
virtual const char * Name() const
Definition: fe_coll.hpp:264
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:992
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:757
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:675
virtual ~L2_FECollection()
Definition: fe_coll.cpp:1281
virtual const char * Name() const
Definition: fe_coll.hpp:490
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:296
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:917
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:520
void InitFaces(const int p, const int dim, const int map_type)
Definition: fe_coll.cpp:1323
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:578
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:556
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:486
int GetOrder() const
Returns the order of the finite element.
Definition: fe.hpp:89
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:221
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1020
void Reset() const
Definition: fe.hpp:1944
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:857
ND_FECollection(const int p, const int dim)
Definition: fe_coll.cpp:1446
Class for quadratic FE on triangle.
Definition: fe.hpp:474
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:400
Class for quadratic FE on interval.
Definition: fe.hpp:445
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:864
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:960
virtual const char * Name() const
Definition: fe_coll.hpp:142
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:139
L2_FECollection(const int p, const int dim, const int type=0)
Definition: fe_coll.cpp:1175
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:709
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:695
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:478
virtual const char * Name() const
Definition: fe_coll.hpp:214
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:896
RT_FECollection(const int p, const int dim, const int map_type)
Definition: fe_coll.hpp:131
virtual const char * Name() const
Definition: fe_coll.hpp:372
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:793
virtual const char * Name() const
Definition: fe_coll.hpp:569
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:425
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:646
Finite element collection on a macro-element.
Definition: fe_coll.hpp:636
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:373
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:134
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:157
Class for refined bi-linear FE on quadrilateral.
Definition: fe.hpp:1029
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:169
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:412
virtual const char * Name() const
Definition: fe_coll.hpp:174
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:548
virtual ~RT_FECollection()
Definition: fe_coll.cpp:1425
virtual int DofForGeometry(int _GeomType) const
Definition: fe_coll.hpp:736
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:464
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1406
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:976
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:288
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:512
virtual const char * Name() const
Definition: fe_coll.hpp:393
Class for bilinear FE on quad with nodes at the 4 Gaussian points.
Definition: fe.hpp:420
Class for refined linear FE on interval.
Definition: fe.hpp:976
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1005
Class for refined linear FE on triangle.
Definition: fe.hpp:996
Class for refined linear FE on tetrahedron.
Definition: fe.hpp:1016
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:123
Class for constant FE on triangle.
Definition: fe.hpp:613
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:223
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:498
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:150
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:565
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:769
Class for refined trilinear FE on a hexahedron.
Definition: fe.hpp:1049
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:724
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:198
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:515
virtual ~ND_FECollection()
Definition: fe_coll.cpp:1590
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:785
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:333
Class for quadratic FE on tetrahedron.
Definition: fe.hpp:668
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:738
RT_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL)
Definition: fe_coll.cpp:1434
H1_FECollection(const int p, const int dim=3, const int type=0)
Definition: fe_coll.cpp:1039
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:598
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:452
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:268
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:717
virtual ~Local_FECollection()
Definition: fe_coll.hpp:742
virtual const char * Name() const
Definition: fe_coll.hpp:476
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:880
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1165
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:118
Class for linear FE on tetrahedron.
Definition: fe.hpp:643
Class for linear FE on interval.
Definition: fe.hpp:342
Crouzeix-Raviart finite element on quadrilateral.
Definition: fe.hpp:716
Class for linear FE on triangle.
Definition: fe.hpp:362
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:946
Class for quadratic FE on triangle with nodes at the "Gaussian" points.
Definition: fe.hpp:498
virtual const char * Name() const
Definition: fe_coll.hpp:68
virtual const char * Name() const
Definition: fe_coll.hpp:611
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:332
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:265
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:312
NURBSFECollection(int Order)
Definition: fe_coll.hpp:193
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:119
L2_FECollection DG_FECollection
Definition: fe_coll.hpp:116
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:323
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:78
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:309
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:386
virtual const char * Name() const
Definition: fe_coll.hpp:328
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:97
virtual const char * Name() const
Definition: fe_coll.hpp:240
virtual const char * Name() const
Definition: fe_coll.hpp:414
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:106
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:805
virtual const char * Name() const
Definition: fe_coll.hpp:350
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:612
virtual const char * Name() const
Definition: fe_coll.hpp:677
int GetDof() const
Returns the degrees of freedom in the FE space.
Definition: fe.hpp:86
Class for bi-quadratic FE on quadrilateral.
Definition: fe.hpp:514
Class for tri-linear FE on cube.
Definition: fe.hpp:681
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:620
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:189
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:533
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:137
virtual const char * Name() const
Definition: fe_coll.hpp:511
Class for bilinear FE on quadrilateral.
Definition: fe.hpp:384
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:841
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:683
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:660
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:459
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1654
virtual int * DofOrderForOrientation(int GeomType, int Or) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:104
static FiniteElementCollection * New(const char *name)
Definition: fe_coll.cpp:38
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:409
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:821
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:905
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:171
virtual const char * Name() const
Definition: fe_coll.hpp:36
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:173
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:63
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
Definition: fe.hpp:553
virtual const char * Name() const
Definition: fe_coll.hpp:434
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:439
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1563
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1146
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:360
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:65
virtual const char * Name() const
Definition: fe_coll.hpp:699
Tensor products of 1D FEs (only degree 2 is functional)
Definition: fe.hpp:955
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:244
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:95
int HasFaceDofs(int GeomType) const
Definition: fe_coll.cpp:25
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:158
void Reset() const
Definition: fe_coll.hpp:200
Class for cubic FE on tetrahedron.
Definition: fe.hpp:600
virtual const char * Name() const
Definition: fe_coll.hpp:740
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:747
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1266
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:52
virtual const char * Name() const
Definition: fe_coll.hpp:720
virtual const char * Name() const
Definition: fe_coll.hpp:527
virtual const FiniteElement * FiniteElementForGeometry(int _GeomType) const
Definition: fe_coll.hpp:734
virtual int DofForGeometry(int GeomType) const =0
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:525
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1667
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:1600
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:631
virtual const char * Name() const
Definition: fe_coll.hpp:284
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:731
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:480
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:249
virtual ~FiniteElementCollection()
Definition: fe_coll.hpp:46
virtual ~NURBSFECollection()
Definition: fe_coll.hpp:216
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:346
Crouzeix-Raviart finite element on triangle.
Definition: fe.hpp:704
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:233
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const =0
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:83
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:586