MFEM  v4.0
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 
22 /** Collection of finite elements from the same family in multiple dimensions.
23  This class is used to match the degrees of freedom of a FiniteElementSpace
24  between elements, and to provide the finite element restriction from an
25  element to its boundary. */
27 {
28 protected:
29  template <Geometry::Type geom>
30  static inline void GetNVE(int &nv, int &ne);
31 
32  template <Geometry::Type geom, typename v_t>
33  static inline void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo,
34  const int edge_info);
35 
36  template <Geometry::Type geom, Geometry::Type f_geom,
37  typename v_t, typename e_t, typename eo_t>
38  static inline void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo,
39  int &nf, int &f, Geometry::Type &fg, int &fo,
40  const int face_info);
41 
42 public:
43  virtual const FiniteElement *
44  FiniteElementForGeometry(Geometry::Type GeomType) const = 0;
45 
46  virtual int DofForGeometry(Geometry::Type GeomType) const = 0;
47 
48  /** @brief Returns an array, say p, that maps a local permuted index i to
49  a local base index: base_i = p[i]. */
50  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
51  int Or) const = 0;
52 
53  virtual const char * Name() const { return "Undefined"; }
54 
55  int HasFaceDofs(Geometry::Type GeomType) const;
56 
58  Geometry::Type GeomType) const
59  {
60  return FiniteElementForGeometry(GeomType);
61  }
62 
64 
66 
67  /** @brief Factory method: return a newly allocated FiniteElementCollection
68  according to the given name. */
69  static FiniteElementCollection *New(const char *name);
70 
71  /** @brief Get the local dofs for a given sub-manifold.
72 
73  Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex,
74  1D - edge, 2D - face) including those on its boundary. The local index of
75  the sub-manifold (inside Geom) and its orientation are given by the
76  parameter Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed
77  that 0 <= SDim <= Dim(Geom). */
78  void SubDofOrder(Geometry::Type Geom, int SDim, int Info,
79  Array<int> &dofs) const;
80 };
81 
82 /// Arbitrary order H1-conforming (continuous) finite elements.
84 {
85 
86 protected:
87  int b_type;
88  char h1_name[32];
91  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
92 
93 public:
94  explicit H1_FECollection(const int p, const int dim = 3,
95  const int btype = BasisType::GaussLobatto);
96 
98  Geometry::Type GeomType) const
99  { return H1_Elements[GeomType]; }
100  virtual int DofForGeometry(Geometry::Type GeomType) const
101  { return H1_dof[GeomType]; }
102  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
103  int Or) const;
104  virtual const char *Name() const { return h1_name; }
106 
107  int GetBasisType() const { return b_type; }
108  /// Get the Cartesian to local H1 dof map
109  const int *GetDofMap(Geometry::Type GeomType) const;
110 
111  virtual ~H1_FECollection();
112 };
113 
114 /** Arbitrary order H1-conforming (continuous) finite elements with positive
115  basis functions. */
117 {
118 public:
119  explicit H1Pos_FECollection(const int p, const int dim = 3)
120  : H1_FECollection(p, dim, BasisType::Positive) { }
121 };
122 
123 /** Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the
124  interface between mesh elements (faces,edges,vertices); these are the trace
125  FEs of the H1-conforming FEs. */
127 {
128 public:
129  H1_Trace_FECollection(const int p, const int dim,
130  const int btype = BasisType::GaussLobatto);
131 };
132 
133 /// Arbitrary order "L2-conforming" discontinuous finite elements.
135 {
136 private:
137  int b_type; // BasisType
138  char d_name[32];
141  int *SegDofOrd[2]; // for rotating segment dofs in 1D
142  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
143  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
144 
145 public:
146  L2_FECollection(const int p, const int dim,
147  const int btype = BasisType::GaussLegendre,
148  const int map_type = FiniteElement::VALUE);
149 
151  Geometry::Type GeomType) const
152  { return L2_Elements[GeomType]; }
153  virtual int DofForGeometry(Geometry::Type GeomType) const
154  {
155  if (L2_Elements[GeomType])
156  {
157  return L2_Elements[GeomType]->GetDof();
158  }
159  return 0;
160  }
161  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
162  int Or) const;
163  virtual const char *Name() const { return d_name; }
164 
166  Geometry::Type GeomType) const
167  {
168  return Tr_Elements[GeomType];
169  }
170 
171  int GetBasisType() const { return b_type; }
172 
173  virtual ~L2_FECollection();
174 };
175 
176 /// Declare an alternative name for L2_FECollection = DG_FECollection
178 
179 /// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
181 {
182 protected:
183  int ob_type; // open BasisType
184  char rt_name[32];
187  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
188 
189  // Initialize only the face elements
190  void InitFaces(const int p, const int dim, const int map_type,
191  const bool signs);
192 
193  // Constructor used by the constructor of the RT_Trace_FECollection and
194  // DG_Interface_FECollection classes
195  RT_FECollection(const int p, const int dim, const int map_type,
196  const bool signs,
197  const int ob_type = BasisType::GaussLegendre);
198 
199 public:
200  RT_FECollection(const int p, const int dim,
201  const int cb_type = BasisType::GaussLobatto,
202  const int ob_type = BasisType::GaussLegendre);
203 
205  Geometry::Type GeomType) const
206  { return RT_Elements[GeomType]; }
207  virtual int DofForGeometry(Geometry::Type GeomType) const
208  { return RT_dof[GeomType]; }
209  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
210  int Or) const;
211  virtual const char *Name() const { return rt_name; }
213 
214  virtual ~RT_FECollection();
215 };
216 
217 /** Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the
218  interface between mesh elements (faces); these are the normal trace FEs of
219  the H(div)-conforming FEs. */
221 {
222 public:
223  RT_Trace_FECollection(const int p, const int dim,
224  const int map_type = FiniteElement::INTEGRAL,
225  const int ob_type = BasisType::GaussLegendre);
226 };
227 
228 /** Arbitrary order discontinuous finite elements defined on the interface
229  between mesh elements (faces). The functions in this space are single-valued
230  on each face and are discontinuous across its boundary. */
232 {
233 public:
234  DG_Interface_FECollection(const int p, const int dim,
235  const int map_type = FiniteElement::VALUE,
236  const int ob_type = BasisType::GaussLegendre);
237 };
238 
239 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
241 {
242 protected:
243  char nd_name[32];
246  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
247 
248 public:
249  ND_FECollection(const int p, const int dim,
250  const int cb_type = BasisType::GaussLobatto,
251  const int ob_type = BasisType::GaussLegendre);
252 
254  const
255  { return ND_Elements[GeomType]; }
256  virtual int DofForGeometry(Geometry::Type GeomType) const
257  { return ND_dof[GeomType]; }
258  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
259  int Or) const;
260  virtual const char *Name() const { return nd_name; }
262 
263  virtual ~ND_FECollection();
264 };
265 
266 /** Arbitrary order H(curl)-trace finite elements defined on the interface
267  between mesh elements (faces,edges); these are the tangential trace FEs of
268  the H(curl)-conforming FEs. */
270 {
271 public:
272  ND_Trace_FECollection(const int p, const int dim,
273  const int cb_type = BasisType::GaussLobatto,
274  const int ob_type = BasisType::GaussLegendre);
275 };
276 
277 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
279 {
280 private:
281  NURBS1DFiniteElement *SegmentFE;
282  NURBS2DFiniteElement *QuadrilateralFE;
283  NURBS3DFiniteElement *ParallelepipedFE;
284 
285  mutable int mOrder; // >= 1 or VariableOrder
286  // The 'name' can be:
287  // 1) name = "NURBS" + "number", for fixed order, or
288  // 2) name = "NURBS", for VariableOrder.
289  // The name is updated before writing it to a stream, for example, see
290  // FiniteElementSpace::Save().
291  mutable char name[16];
292 
293 public:
294  enum { VariableOrder = -1 };
295 
296  /** @brief The parameter @a Order must be either a positive number, for fixed
297  order, or VariableOrder (default). */
298  explicit NURBSFECollection(int Order = VariableOrder);
299 
300  void Reset() const
301  {
302  SegmentFE->Reset();
303  QuadrilateralFE->Reset();
304  ParallelepipedFE->Reset();
305  }
306 
307  /** @brief Get the order of the NURBS collection: either a positive number,
308  when using fixed order, or VariableOrder. */
309  int GetOrder() const { return mOrder; }
310 
311  /** @brief Set the order and the name, based on the given @a Order: either a
312  positive number for fixed order, or VariableOrder. */
313  void SetOrder(int Order) const;
314 
315  virtual const FiniteElement *
317 
318  virtual int DofForGeometry(Geometry::Type GeomType) const;
319 
320  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
321  int Or) const;
322 
323  virtual const char *Name() const { return name; }
324 
326 
327  virtual ~NURBSFECollection();
328 };
329 
330 
331 /// Piecewise-(bi)linear continuous finite elements.
333 {
334 private:
335  const PointFiniteElement PointFE;
336  const Linear1DFiniteElement SegmentFE;
337  const Linear2DFiniteElement TriangleFE;
338  const BiLinear2DFiniteElement QuadrilateralFE;
339  const Linear3DFiniteElement TetrahedronFE;
340  const TriLinear3DFiniteElement ParallelepipedFE;
341  const H1_WedgeElement WedgeFE;
342 public:
343  LinearFECollection() : WedgeFE(1) { }
344 
345  virtual const FiniteElement *
347 
348  virtual int DofForGeometry(Geometry::Type GeomType) const;
349 
350  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
351  int Or) const;
352 
353  virtual const char * Name() const { return "Linear"; }
354 };
355 
356 /// Piecewise-(bi)quadratic continuous finite elements.
358 {
359 private:
360  const PointFiniteElement PointFE;
361  const Quad1DFiniteElement SegmentFE;
362  const Quad2DFiniteElement TriangleFE;
363  const BiQuad2DFiniteElement QuadrilateralFE;
364  const Quadratic3DFiniteElement TetrahedronFE;
365  const LagrangeHexFiniteElement ParallelepipedFE;
366  const H1_WedgeElement WedgeFE;
367 
368 public:
369  QuadraticFECollection() : ParallelepipedFE(2), WedgeFE(2) { }
370 
371  virtual const FiniteElement *
373 
374  virtual int DofForGeometry(Geometry::Type GeomType) const;
375 
376  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
377  int Or) const;
378 
379  virtual const char * Name() const { return "Quadratic"; }
380 };
381 
382 /// Version of QuadraticFECollection with positive basis functions.
384 {
385 private:
386  const QuadPos1DFiniteElement SegmentFE;
387  const BiQuadPos2DFiniteElement QuadrilateralFE;
388 
389 public:
391 
392  virtual const FiniteElement *
394 
395  virtual int DofForGeometry(Geometry::Type GeomType) const;
396 
397  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
398  int Or) const;
399 
400  virtual const char * Name() const { return "QuadraticPos"; }
401 };
402 
403 /// Piecewise-(bi)cubic continuous finite elements.
405 {
406 private:
407  const PointFiniteElement PointFE;
408  const Cubic1DFiniteElement SegmentFE;
409  const Cubic2DFiniteElement TriangleFE;
410  const BiCubic2DFiniteElement QuadrilateralFE;
411  const Cubic3DFiniteElement TetrahedronFE;
412  const LagrangeHexFiniteElement ParallelepipedFE;
413  const H1_WedgeElement WedgeFE;
414 
415 public:
417  : ParallelepipedFE(3), WedgeFE(3, BasisType::ClosedUniform) { }
418 
419  virtual const FiniteElement *
421 
422  virtual int DofForGeometry(Geometry::Type GeomType) const;
423 
424  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
425  int Or) const;
426 
427  virtual const char * Name() const { return "Cubic"; }
428 };
429 
430 /// Crouzeix-Raviart nonconforming elements in 2D.
432 {
433 private:
434  const P0SegmentFiniteElement SegmentFE;
435  const CrouzeixRaviartFiniteElement TriangleFE;
436  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
437 public:
438  CrouzeixRaviartFECollection() : SegmentFE(1) { }
439 
440  virtual const FiniteElement *
442 
443  virtual int DofForGeometry(Geometry::Type GeomType) const;
444 
445  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
446  int Or) const;
447 
448  virtual const char * Name() const { return "CrouzeixRaviart"; }
449 };
450 
451 /// Piecewise-linear nonconforming finite elements in 3D.
453 {
454 private:
455  const P0TriangleFiniteElement TriangleFE;
456  const P1TetNonConfFiniteElement TetrahedronFE;
457  const P0QuadFiniteElement QuadrilateralFE;
458  const RotTriLinearHexFiniteElement ParallelepipedFE;
459 
460 public:
462 
463  virtual const FiniteElement *
465 
466  virtual int DofForGeometry(Geometry::Type GeomType) const;
467 
468  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
469  int Or) const;
470 
471  virtual const char * Name() const { return "LinearNonConf3D"; }
472 };
473 
474 
475 /** First order Raviart-Thomas finite elements in 2D. This class is kept only
476  for backward compatibility, consider using RT_FECollection instead. */
478 {
479 private:
480  const P0SegmentFiniteElement SegmentFE; // normal component on edge
481  const RT0TriangleFiniteElement TriangleFE;
482  const RT0QuadFiniteElement QuadrilateralFE;
483 public:
484  RT0_2DFECollection() : SegmentFE(0) { }
485 
486  virtual const FiniteElement *
488 
489  virtual int DofForGeometry(Geometry::Type GeomType) const;
490 
491  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
492  int Or) const;
493 
494  virtual const char * Name() const { return "RT0_2D"; }
495 };
496 
497 /** Second order Raviart-Thomas finite elements in 2D. This class is kept only
498  for backward compatibility, consider using RT_FECollection instead. */
500 {
501 private:
502  const P1SegmentFiniteElement SegmentFE; // normal component on edge
503  const RT1TriangleFiniteElement TriangleFE;
504  const RT1QuadFiniteElement QuadrilateralFE;
505 public:
507 
508  virtual const FiniteElement *
510 
511  virtual int DofForGeometry(Geometry::Type GeomType) const;
512 
513  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
514  int Or) const;
515 
516  virtual const char * Name() const { return "RT1_2D"; }
517 };
518 
519 /** Third order Raviart-Thomas finite elements in 2D. This class is kept only
520  for backward compatibility, consider using RT_FECollection instead. */
522 {
523 private:
524  const P2SegmentFiniteElement SegmentFE; // normal component on edge
525  const RT2TriangleFiniteElement TriangleFE;
526  const RT2QuadFiniteElement QuadrilateralFE;
527 public:
529 
530  virtual const FiniteElement *
532 
533  virtual int DofForGeometry(Geometry::Type GeomType) const;
534 
535  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
536  int Or) const;
537 
538  virtual const char * Name() const { return "RT2_2D"; }
539 };
540 
541 /** Piecewise-constant discontinuous finite elements in 2D. This class is kept
542  only for backward compatibility, consider using L2_FECollection instead. */
544 {
545 private:
546  const P0TriangleFiniteElement TriangleFE;
547  const P0QuadFiniteElement QuadrilateralFE;
548 public:
550 
551  virtual const FiniteElement *
553 
554  virtual int DofForGeometry(Geometry::Type GeomType) const;
555 
556  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
557  int Or) const;
558 
559  virtual const char * Name() const { return "Const2D"; }
560 };
561 
562 /** Piecewise-linear discontinuous finite elements in 2D. This class is kept
563  only for backward compatibility, consider using L2_FECollection instead. */
565 {
566 private:
567  const Linear2DFiniteElement TriangleFE;
568  const BiLinear2DFiniteElement QuadrilateralFE;
569 
570 public:
572 
573  virtual const FiniteElement *
575 
576  virtual int DofForGeometry(Geometry::Type GeomType) const;
577 
578  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
579  int Or) const;
580 
581  virtual const char * Name() const { return "LinearDiscont2D"; }
582 };
583 
584 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
586 {
587 private:
588  // const CrouzeixRaviartFiniteElement TriangleFE;
589  const GaussLinear2DFiniteElement TriangleFE;
590  const GaussBiLinear2DFiniteElement QuadrilateralFE;
591 
592 public:
594 
595  virtual const FiniteElement *
597 
598  virtual int DofForGeometry(Geometry::Type GeomType) const;
599 
600  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
601  int Or) const;
602 
603  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
604 };
605 
606 /// Linear (P1) finite elements on quadrilaterals.
608 {
609 private:
610  const P1OnQuadFiniteElement QuadrilateralFE;
611 public:
613  virtual const FiniteElement *
615  virtual int DofForGeometry(Geometry::Type GeomType) const;
616  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
617  int Or) const;
618  virtual const char * Name() const { return "P1OnQuad"; }
619 };
620 
621 /** Piecewise-quadratic discontinuous finite elements in 2D. This class is kept
622  only for backward compatibility, consider using L2_FECollection instead. */
624 {
625 private:
626  const Quad2DFiniteElement TriangleFE;
627  const BiQuad2DFiniteElement QuadrilateralFE;
628 
629 public:
631 
632  virtual const FiniteElement *
634 
635  virtual int DofForGeometry(Geometry::Type GeomType) const;
636 
637  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
638  int Or) const;
639 
640  virtual const char * Name() const { return "QuadraticDiscont2D"; }
641 };
642 
643 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
645 {
646 private:
647  const BiQuadPos2DFiniteElement QuadrilateralFE;
648 
649 public:
651  virtual const FiniteElement *
653  virtual int DofForGeometry(Geometry::Type GeomType) const;
654  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
655  int Or) const
656  { return NULL; }
657  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
658 };
659 
660 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
662 {
663 private:
664  // const Quad2DFiniteElement TriangleFE;
665  const GaussQuad2DFiniteElement TriangleFE;
666  const GaussBiQuad2DFiniteElement QuadrilateralFE;
667 
668 public:
670 
671  virtual const FiniteElement *
673 
674  virtual int DofForGeometry(Geometry::Type GeomType) const;
675 
676  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
677  int Or) const;
678 
679  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
680 };
681 
682 /** Piecewise-cubic discontinuous finite elements in 2D. This class is kept
683  only for backward compatibility, consider using L2_FECollection instead. */
685 {
686 private:
687  const Cubic2DFiniteElement TriangleFE;
688  const BiCubic2DFiniteElement QuadrilateralFE;
689 
690 public:
692 
693  virtual const FiniteElement *
695 
696  virtual int DofForGeometry(Geometry::Type GeomType) const;
697 
698  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
699  int Or) const;
700 
701  virtual const char * Name() const { return "CubicDiscont2D"; }
702 };
703 
704 /** Piecewise-constant discontinuous finite elements in 3D. This class is kept
705  only for backward compatibility, consider using L2_FECollection instead. */
707 {
708 private:
709  const P0TetFiniteElement TetrahedronFE;
710  const P0HexFiniteElement ParallelepipedFE;
711  const L2_WedgeElement WedgeFE;
712 
713 public:
714  Const3DFECollection() : WedgeFE(0) { }
715 
716  virtual const FiniteElement *
718 
719  virtual int DofForGeometry(Geometry::Type GeomType) const;
720 
721  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
722  int Or) const;
723 
724  virtual const char * Name() const { return "Const3D"; }
725 };
726 
727 /** Piecewise-linear discontinuous finite elements in 3D. This class is kept
728  only for backward compatibility, consider using L2_FECollection instead. */
730 {
731 private:
732  const Linear3DFiniteElement TetrahedronFE;
733  const TriLinear3DFiniteElement ParallelepipedFE;
734 
735 public:
737 
738  virtual const FiniteElement *
740 
741  virtual int DofForGeometry(Geometry::Type GeomType) const;
742 
743  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
744  int Or) const;
745 
746  virtual const char * Name() const { return "LinearDiscont3D"; }
747 };
748 
749 /** Piecewise-quadratic discontinuous finite elements in 3D. This class is kept
750  only for backward compatibility, consider using L2_FECollection instead. */
752 {
753 private:
754  const Quadratic3DFiniteElement TetrahedronFE;
755  const LagrangeHexFiniteElement ParallelepipedFE;
756 
757 public:
758  QuadraticDiscont3DFECollection () : ParallelepipedFE(2) { }
759 
760  virtual const FiniteElement *
762 
763  virtual int DofForGeometry(Geometry::Type GeomType) const;
764 
765  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
766  int Or) const;
767 
768  virtual const char * Name() const { return "QuadraticDiscont3D"; }
769 };
770 
771 /// Finite element collection on a macro-element.
773 {
774 private:
775  const PointFiniteElement PointFE;
776  const RefinedLinear1DFiniteElement SegmentFE;
777  const RefinedLinear2DFiniteElement TriangleFE;
778  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
779  const RefinedLinear3DFiniteElement TetrahedronFE;
780  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
781 
782 public:
784 
785  virtual const FiniteElement *
787 
788  virtual int DofForGeometry(Geometry::Type GeomType) const;
789 
790  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
791  int Or) const;
792 
793  virtual const char * Name() const { return "RefinedLinear"; }
794 };
795 
796 /** Lowest order Nedelec finite elements in 3D. This class is kept only for
797  backward compatibility, consider using the new ND_FECollection instead. */
799 {
800 private:
801  const Nedelec1HexFiniteElement HexahedronFE;
802  const Nedelec1TetFiniteElement TetrahedronFE;
803 
804 public:
806 
807  virtual const FiniteElement *
809 
810  virtual int DofForGeometry(Geometry::Type GeomType) const;
811 
812  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
813  int Or) const;
814 
815  virtual const char * Name() const { return "ND1_3D"; }
816 };
817 
818 /** First order Raviart-Thomas finite elements in 3D. This class is kept only
819  for backward compatibility, consider using RT_FECollection instead. */
821 {
822 private:
823  const P0TriangleFiniteElement TriangleFE;
824  const P0QuadFiniteElement QuadrilateralFE;
825  const RT0HexFiniteElement HexahedronFE;
826  const RT0TetFiniteElement TetrahedronFE;
827 public:
829 
830  virtual const FiniteElement *
832 
833  virtual int DofForGeometry(Geometry::Type GeomType) const;
834 
835  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
836  int Or) const;
837 
838  virtual const char * Name() const { return "RT0_3D"; }
839 };
840 
841 /** Second order Raviart-Thomas finite elements in 3D. This class is kept only
842  for backward compatibility, consider using RT_FECollection instead. */
844 {
845 private:
846  const Linear2DFiniteElement TriangleFE;
847  const BiLinear2DFiniteElement QuadrilateralFE;
848  const RT1HexFiniteElement HexahedronFE;
849 public:
851 
852  virtual const FiniteElement *
854 
855  virtual int DofForGeometry(Geometry::Type GeomType) const;
856 
857  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
858  int Or) const;
859 
860  virtual const char * Name() const { return "RT1_3D"; }
861 };
862 
863 /// Discontinuous collection defined locally by a given finite element.
865 {
866 private:
867  char d_name[32];
868  Geometry::Type GeomType;
869  FiniteElement *Local_Element;
870 
871 public:
872  Local_FECollection(const char *fe_name);
873 
875  Geometry::Type _GeomType) const
876  { return (GeomType == _GeomType) ? Local_Element : NULL; }
877  virtual int DofForGeometry(Geometry::Type _GeomType) const
878  { return (GeomType == _GeomType) ? Local_Element->GetDof() : 0; }
879  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
880  int Or) const
881  { return NULL; }
882  virtual const char *Name() const { return d_name; }
883 
884  virtual ~Local_FECollection() { delete Local_Element; }
885 };
886 
887 }
888 
889 #endif
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:834
Abstract class for Finite Elements.
Definition: fe.hpp:229
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:309
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:951
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:278
virtual const char * Name() const
Definition: fe_coll.hpp:427
virtual const char * Name() const
Definition: fe_coll.hpp:768
virtual const char * Name() const
Definition: fe_coll.hpp:679
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1728
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1165
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1157
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:661
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:186
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1027
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:973
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:244
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.hpp:654
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1206
Linear 1D element with nodes 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:1350
Quadratic 1D element with nodes the Gaussian points in [0,1] (trace of RT2)
Definition: fe.hpp:1360
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1122
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1236
virtual const char * Name() const
Definition: fe_coll.hpp:793
For scalar fields; preserves volume integrals.
Definition: fe.hpp:274
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:724
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:97
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:2505
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1899
RT_FECollection(const int p, const int dim, const int map_type, const bool signs, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:1990
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:570
virtual const char * Name() const
Definition: fe_coll.hpp:581
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2121
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2539
static const int NumGeom
Definition: geom.hpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:379
const Geometry::Type geom
Definition: ex1.cpp:40
virtual ~L2_FECollection()
Definition: fe_coll.cpp:1915
virtual const char * Name() const
Definition: fe_coll.hpp:618
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2513
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:755
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const =0
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:612
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:332
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1658
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:924
void Reset() const
Definition: fe.hpp:2865
Class for quadratic FE on triangle.
Definition: fe.hpp:933
Class for quadratic FE on interval.
Definition: fe.hpp:904
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1373
virtual const char * Name() const
Definition: fe_coll.hpp:211
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2389
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:989
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:658
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:100
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:2003
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2139
virtual const char * Name() const
Definition: fe_coll.hpp:323
DG_Interface_FECollection(const int p, const int dim, const int map_type=FiniteElement::VALUE, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2187
int HasFaceDofs(Geometry::Type GeomType) const
Definition: fe_coll.cpp:25
virtual const char * Name() const
Definition: fe_coll.hpp:494
virtual const char * Name() const
Definition: fe_coll.hpp:701
Finite element collection on a macro-element.
Definition: fe_coll.hpp:772
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:916
Class for refined bi-linear FE on quadrilateral.
Definition: fe.hpp:1489
Possible basis types. Note that not all elements can use all BasisType(s).
Definition: fe.hpp:27
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default)...
Definition: fe_coll.cpp:2482
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1180
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:260
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1244
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:622
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:648
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:692
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2155
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1748
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:150
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:404
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:820
virtual const char * Name() const
Definition: fe_coll.hpp:516
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1695
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:875
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:536
int dim
Definition: ex3.cpp:48
Class for bilinear FE on quad with nodes at the 4 Gaussian points.
Definition: fe.hpp:879
Class for refined linear FE on interval.
Definition: fe.hpp:1436
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:90
Class for refined linear FE on triangle.
Definition: fe.hpp:1456
Class for refined linear FE on tetrahedron.
Definition: fe.hpp:1476
static void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo, int &nf, int &f, Geometry::Type &fg, int &fo, const int face_info)
Definition: fe_coll.cpp:319
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:185
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:902
Class for constant FE on triangle.
Definition: fe.hpp:1072
virtual int DofForGeometry(Geometry::Type _GeomType) const
Definition: fe_coll.hpp:877
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:289
int GetBasisType() const
Definition: fe_coll.hpp:171
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1361
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1351
int GetBasisType() const
Definition: fe_coll.hpp:107
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1502
ND_Trace_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2421
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2207
Class for refined trilinear FE on a hexahedron.
Definition: fe.hpp:1509
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:864
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:644
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2409
Class for quadratic FE on tetrahedron.
Definition: fe.hpp:1127
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:807
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:959
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:778
void SubDofOrder(Geometry::Type Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
Definition: fe_coll.cpp:364
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1055
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1257
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:383
virtual ~Local_FECollection()
Definition: fe_coll.hpp:884
virtual const char * Name() const
Definition: fe_coll.hpp:603
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:57
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:578
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1716
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1311
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1273
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.hpp:879
Class for linear FE on tetrahedron.
Definition: fe.hpp:1102
Class for linear FE on interval.
Definition: fe.hpp:802
Crouzeix-Raviart finite element on quadrilateral.
Definition: fe.hpp:1176
Class for linear FE on triangle.
Definition: fe.hpp:822
Class for quadratic FE on triangle with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:957
virtual const char * Name() const
Definition: fe_coll.hpp:104
virtual const char * Name() const
Definition: fe_coll.hpp:746
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1676
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:207
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:89
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:452
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:431
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:765
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1294
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1106
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:180
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:997
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:177
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:119
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:153
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:299
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:165
virtual const char * Name() const
Definition: fe_coll.hpp:448
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1434
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const =0
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const char * Name() const
Definition: fe_coll.hpp:353
virtual const char * Name() const
Definition: fe_coll.hpp:538
void SetOrder(int Order) const
Set the order and the name, based on the given Order: either a positive number for fixed order...
Definition: fe_coll.cpp:2492
virtual const char * Name() const
Definition: fe_coll.hpp:471
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:861
virtual const char * Name() const
Definition: fe_coll.hpp:815
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:311
Class for bi-quadratic FE on quadrilateral.
Definition: fe.hpp:973
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2366
Class for tri-linear FE on cube.
Definition: fe.hpp:1140
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type _GeomType) const
Definition: fe_coll.hpp:874
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:741
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:595
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1281
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1404
virtual const char * Name() const
Definition: fe_coll.hpp:640
Class for bilinear FE on quadrilateral.
Definition: fe.hpp:844
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1335
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:585
virtual const char * Name() const
Definition: fe_coll.hpp:163
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:47
Class for linear FE on triangle with nodes at the 3 &quot;Gaussian&quot; points.
Definition: fe.hpp:868
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1318
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1019
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:727
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1006
virtual const char * Name() const
Definition: fe_coll.hpp:53
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1465
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:256
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1480
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
Definition: fe.hpp:1012
virtual const char * Name() const
Definition: fe_coll.hpp:559
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2532
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:635
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1063
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1075
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:204
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1219
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:792
virtual const char * Name() const
Definition: fe_coll.hpp:838
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:675
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1418
RT_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2167
Tensor products of 1D FEs (only degree 2 is functional)
Definition: fe.hpp:1415
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:553
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:357
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1196
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:848
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:240
void Reset() const
Definition: fe_coll.hpp:300
Class for cubic FE on tetrahedron.
Definition: fe.hpp:1059
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1040
For scalar fields; preserves point values.
Definition: fe.hpp:273
virtual const char * Name() const
Definition: fe_coll.hpp:882
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
virtual const char * Name() const
Definition: fe_coll.hpp:860
virtual const char * Name() const
Definition: fe_coll.hpp:657
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1091
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2526
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:937
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:890
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:2440
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1143
virtual const char * Name() const
Definition: fe_coll.hpp:400
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:607
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1389
virtual ~FiniteElementCollection()
Definition: fe_coll.hpp:65
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1130
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:245
Crouzeix-Raviart finite element on triangle.
Definition: fe.hpp:1164
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1452
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:134
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:253