MFEM  v4.1.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-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
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 H1-conforming (continuous) serendipity finite elements;
124  Current implementation works in 2D only; 3D version is in development. */
126 {
127 public:
128  explicit H1Ser_FECollection(const int p, const int dim = 2)
129  : H1_FECollection(p, dim, BasisType::Serendipity) { };
130 };
131 
132 /** Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the
133  interface between mesh elements (faces,edges,vertices); these are the trace
134  FEs of the H1-conforming FEs. */
136 {
137 public:
138  H1_Trace_FECollection(const int p, const int dim,
139  const int btype = BasisType::GaussLobatto);
140 };
141 
142 /// Arbitrary order "L2-conforming" discontinuous finite elements.
144 {
145 private:
146  int b_type; // BasisType
147  char d_name[32];
150  int *SegDofOrd[2]; // for rotating segment dofs in 1D
151  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
152  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
153 
154 public:
155  L2_FECollection(const int p, const int dim,
156  const int btype = BasisType::GaussLegendre,
157  const int map_type = FiniteElement::VALUE);
158 
160  Geometry::Type GeomType) const
161  {
162  return L2_Elements[GeomType];
163  }
164  virtual int DofForGeometry(Geometry::Type GeomType) const
165  {
166  if (L2_Elements[GeomType])
167  {
168  return L2_Elements[GeomType]->GetDof();
169  }
170  return 0;
171  }
172  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
173  int Or) const;
174  virtual const char *Name() const { return d_name; }
175 
177  Geometry::Type GeomType) const
178  {
179  return Tr_Elements[GeomType];
180  }
181 
182  int GetBasisType() const { return b_type; }
183 
184  virtual ~L2_FECollection();
185 };
186 
187 /// Declare an alternative name for L2_FECollection = DG_FECollection
189 
190 /// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
192 {
193 protected:
194  int ob_type; // open BasisType
195  char rt_name[32];
198  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
199 
200  // Initialize only the face elements
201  void InitFaces(const int p, const int dim, const int map_type,
202  const bool signs);
203 
204  // Constructor used by the constructor of the RT_Trace_FECollection and
205  // DG_Interface_FECollection classes
206  RT_FECollection(const int p, const int dim, const int map_type,
207  const bool signs,
208  const int ob_type = BasisType::GaussLegendre);
209 
210 public:
211  RT_FECollection(const int p, const int dim,
212  const int cb_type = BasisType::GaussLobatto,
213  const int ob_type = BasisType::GaussLegendre);
214 
216  Geometry::Type GeomType) const
217  { return RT_Elements[GeomType]; }
218  virtual int DofForGeometry(Geometry::Type GeomType) const
219  { return RT_dof[GeomType]; }
220  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
221  int Or) const;
222  virtual const char *Name() const { return rt_name; }
224 
225  virtual ~RT_FECollection();
226 };
227 
228 /** Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the
229  interface between mesh elements (faces); these are the normal trace FEs of
230  the H(div)-conforming FEs. */
232 {
233 public:
234  RT_Trace_FECollection(const int p, const int dim,
235  const int map_type = FiniteElement::INTEGRAL,
236  const int ob_type = BasisType::GaussLegendre);
237 };
238 
239 /** Arbitrary order discontinuous finite elements defined on the interface
240  between mesh elements (faces). The functions in this space are single-valued
241  on each face and are discontinuous across its boundary. */
243 {
244 public:
245  DG_Interface_FECollection(const int p, const int dim,
246  const int map_type = FiniteElement::VALUE,
247  const int ob_type = BasisType::GaussLegendre);
248 };
249 
250 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
252 {
253 protected:
254  char nd_name[32];
257  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
258 
259 public:
260  ND_FECollection(const int p, const int dim,
261  const int cb_type = BasisType::GaussLobatto,
262  const int ob_type = BasisType::GaussLegendre);
263 
265  const
266  { return ND_Elements[GeomType]; }
267  virtual int DofForGeometry(Geometry::Type GeomType) const
268  { return ND_dof[GeomType]; }
269  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
270  int Or) const;
271  virtual const char *Name() const { return nd_name; }
273 
274  virtual ~ND_FECollection();
275 };
276 
277 /** Arbitrary order H(curl)-trace finite elements defined on the interface
278  between mesh elements (faces,edges); these are the tangential trace FEs of
279  the H(curl)-conforming FEs. */
281 {
282 public:
283  ND_Trace_FECollection(const int p, const int dim,
284  const int cb_type = BasisType::GaussLobatto,
285  const int ob_type = BasisType::GaussLegendre);
286 };
287 
288 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
290 {
291 private:
292  NURBS1DFiniteElement *SegmentFE;
293  NURBS2DFiniteElement *QuadrilateralFE;
294  NURBS3DFiniteElement *ParallelepipedFE;
295 
296  mutable int mOrder; // >= 1 or VariableOrder
297  // The 'name' can be:
298  // 1) name = "NURBS" + "number", for fixed order, or
299  // 2) name = "NURBS", for VariableOrder.
300  // The name is updated before writing it to a stream, for example, see
301  // FiniteElementSpace::Save().
302  mutable char name[16];
303 
304 public:
305  enum { VariableOrder = -1 };
306 
307  /** @brief The parameter @a Order must be either a positive number, for fixed
308  order, or VariableOrder (default). */
309  explicit NURBSFECollection(int Order = VariableOrder);
310 
311  void Reset() const
312  {
313  SegmentFE->Reset();
314  QuadrilateralFE->Reset();
315  ParallelepipedFE->Reset();
316  }
317 
318  /** @brief Get the order of the NURBS collection: either a positive number,
319  when using fixed order, or VariableOrder. */
320  int GetOrder() const { return mOrder; }
321 
322  /** @brief Set the order and the name, based on the given @a Order: either a
323  positive number for fixed order, or VariableOrder. */
324  void SetOrder(int Order) const;
325 
326  virtual const FiniteElement *
328 
329  virtual int DofForGeometry(Geometry::Type GeomType) const;
330 
331  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
332  int Or) const;
333 
334  virtual const char *Name() const { return name; }
335 
337 
338  virtual ~NURBSFECollection();
339 };
340 
341 
342 /// Piecewise-(bi)linear continuous finite elements.
344 {
345 private:
346  const PointFiniteElement PointFE;
347  const Linear1DFiniteElement SegmentFE;
348  const Linear2DFiniteElement TriangleFE;
349  const BiLinear2DFiniteElement QuadrilateralFE;
350  const Linear3DFiniteElement TetrahedronFE;
351  const TriLinear3DFiniteElement ParallelepipedFE;
352  const H1_WedgeElement WedgeFE;
353 public:
354  LinearFECollection() : WedgeFE(1) { }
355 
356  virtual const FiniteElement *
358 
359  virtual int DofForGeometry(Geometry::Type GeomType) const;
360 
361  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
362  int Or) const;
363 
364  virtual const char * Name() const { return "Linear"; }
365 };
366 
367 /// Piecewise-(bi)quadratic continuous finite elements.
369 {
370 private:
371  const PointFiniteElement PointFE;
372  const Quad1DFiniteElement SegmentFE;
373  const Quad2DFiniteElement TriangleFE;
374  const BiQuad2DFiniteElement QuadrilateralFE;
375  const Quadratic3DFiniteElement TetrahedronFE;
376  const LagrangeHexFiniteElement ParallelepipedFE;
377  const H1_WedgeElement WedgeFE;
378 
379 public:
380  QuadraticFECollection() : ParallelepipedFE(2), WedgeFE(2) { }
381 
382  virtual const FiniteElement *
384 
385  virtual int DofForGeometry(Geometry::Type GeomType) const;
386 
387  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
388  int Or) const;
389 
390  virtual const char * Name() const { return "Quadratic"; }
391 };
392 
393 /// Version of QuadraticFECollection with positive basis functions.
395 {
396 private:
397  const QuadPos1DFiniteElement SegmentFE;
398  const BiQuadPos2DFiniteElement QuadrilateralFE;
399 
400 public:
402 
403  virtual const FiniteElement *
405 
406  virtual int DofForGeometry(Geometry::Type GeomType) const;
407 
408  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
409  int Or) const;
410 
411  virtual const char * Name() const { return "QuadraticPos"; }
412 };
413 
414 /// Piecewise-(bi)cubic continuous finite elements.
416 {
417 private:
418  const PointFiniteElement PointFE;
419  const Cubic1DFiniteElement SegmentFE;
420  const Cubic2DFiniteElement TriangleFE;
421  const BiCubic2DFiniteElement QuadrilateralFE;
422  const Cubic3DFiniteElement TetrahedronFE;
423  const LagrangeHexFiniteElement ParallelepipedFE;
424  const H1_WedgeElement WedgeFE;
425 
426 public:
428  : ParallelepipedFE(3), WedgeFE(3, BasisType::ClosedUniform) { }
429 
430  virtual const FiniteElement *
432 
433  virtual int DofForGeometry(Geometry::Type GeomType) const;
434 
435  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
436  int Or) const;
437 
438  virtual const char * Name() const { return "Cubic"; }
439 };
440 
441 /// Crouzeix-Raviart nonconforming elements in 2D.
443 {
444 private:
445  const P0SegmentFiniteElement SegmentFE;
446  const CrouzeixRaviartFiniteElement TriangleFE;
447  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
448 public:
449  CrouzeixRaviartFECollection() : SegmentFE(1) { }
450 
451  virtual const FiniteElement *
453 
454  virtual int DofForGeometry(Geometry::Type GeomType) const;
455 
456  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
457  int Or) const;
458 
459  virtual const char * Name() const { return "CrouzeixRaviart"; }
460 };
461 
462 /// Piecewise-linear nonconforming finite elements in 3D.
464 {
465 private:
466  const P0TriangleFiniteElement TriangleFE;
467  const P1TetNonConfFiniteElement TetrahedronFE;
468  const P0QuadFiniteElement QuadrilateralFE;
469  const RotTriLinearHexFiniteElement ParallelepipedFE;
470 
471 public:
473 
474  virtual const FiniteElement *
476 
477  virtual int DofForGeometry(Geometry::Type GeomType) const;
478 
479  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
480  int Or) const;
481 
482  virtual const char * Name() const { return "LinearNonConf3D"; }
483 };
484 
485 
486 /** First order Raviart-Thomas finite elements in 2D. This class is kept only
487  for backward compatibility, consider using RT_FECollection instead. */
489 {
490 private:
491  const P0SegmentFiniteElement SegmentFE; // normal component on edge
492  const RT0TriangleFiniteElement TriangleFE;
493  const RT0QuadFiniteElement QuadrilateralFE;
494 public:
495  RT0_2DFECollection() : SegmentFE(0) { }
496 
497  virtual const FiniteElement *
499 
500  virtual int DofForGeometry(Geometry::Type GeomType) const;
501 
502  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
503  int Or) const;
504 
505  virtual const char * Name() const { return "RT0_2D"; }
506 };
507 
508 /** Second order Raviart-Thomas finite elements in 2D. This class is kept only
509  for backward compatibility, consider using RT_FECollection instead. */
511 {
512 private:
513  const P1SegmentFiniteElement SegmentFE; // normal component on edge
514  const RT1TriangleFiniteElement TriangleFE;
515  const RT1QuadFiniteElement QuadrilateralFE;
516 public:
518 
519  virtual const FiniteElement *
521 
522  virtual int DofForGeometry(Geometry::Type GeomType) const;
523 
524  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
525  int Or) const;
526 
527  virtual const char * Name() const { return "RT1_2D"; }
528 };
529 
530 /** Third order Raviart-Thomas finite elements in 2D. This class is kept only
531  for backward compatibility, consider using RT_FECollection instead. */
533 {
534 private:
535  const P2SegmentFiniteElement SegmentFE; // normal component on edge
536  const RT2TriangleFiniteElement TriangleFE;
537  const RT2QuadFiniteElement QuadrilateralFE;
538 public:
540 
541  virtual const FiniteElement *
543 
544  virtual int DofForGeometry(Geometry::Type GeomType) const;
545 
546  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
547  int Or) const;
548 
549  virtual const char * Name() const { return "RT2_2D"; }
550 };
551 
552 /** Piecewise-constant discontinuous finite elements in 2D. This class is kept
553  only for backward compatibility, consider using L2_FECollection instead. */
555 {
556 private:
557  const P0TriangleFiniteElement TriangleFE;
558  const P0QuadFiniteElement QuadrilateralFE;
559 public:
561 
562  virtual const FiniteElement *
564 
565  virtual int DofForGeometry(Geometry::Type GeomType) const;
566 
567  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
568  int Or) const;
569 
570  virtual const char * Name() const { return "Const2D"; }
571 };
572 
573 /** Piecewise-linear discontinuous finite elements in 2D. This class is kept
574  only for backward compatibility, consider using L2_FECollection instead. */
576 {
577 private:
578  const Linear2DFiniteElement TriangleFE;
579  const BiLinear2DFiniteElement QuadrilateralFE;
580 
581 public:
583 
584  virtual const FiniteElement *
586 
587  virtual int DofForGeometry(Geometry::Type GeomType) const;
588 
589  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
590  int Or) const;
591 
592  virtual const char * Name() const { return "LinearDiscont2D"; }
593 };
594 
595 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
597 {
598 private:
599  // const CrouzeixRaviartFiniteElement TriangleFE;
600  const GaussLinear2DFiniteElement TriangleFE;
601  const GaussBiLinear2DFiniteElement QuadrilateralFE;
602 
603 public:
605 
606  virtual const FiniteElement *
608 
609  virtual int DofForGeometry(Geometry::Type GeomType) const;
610 
611  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
612  int Or) const;
613 
614  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
615 };
616 
617 /// Linear (P1) finite elements on quadrilaterals.
619 {
620 private:
621  const P1OnQuadFiniteElement QuadrilateralFE;
622 public:
624  virtual const FiniteElement *
626  virtual int DofForGeometry(Geometry::Type GeomType) const;
627  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
628  int Or) const;
629  virtual const char * Name() const { return "P1OnQuad"; }
630 };
631 
632 /** Piecewise-quadratic discontinuous finite elements in 2D. This class is kept
633  only for backward compatibility, consider using L2_FECollection instead. */
635 {
636 private:
637  const Quad2DFiniteElement TriangleFE;
638  const BiQuad2DFiniteElement QuadrilateralFE;
639 
640 public:
642 
643  virtual const FiniteElement *
645 
646  virtual int DofForGeometry(Geometry::Type GeomType) const;
647 
648  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
649  int Or) const;
650 
651  virtual const char * Name() const { return "QuadraticDiscont2D"; }
652 };
653 
654 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
656 {
657 private:
658  const BiQuadPos2DFiniteElement QuadrilateralFE;
659 
660 public:
662  virtual const FiniteElement *
664  virtual int DofForGeometry(Geometry::Type GeomType) const;
665  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
666  int Or) const
667  { return NULL; }
668  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
669 };
670 
671 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
673 {
674 private:
675  // const Quad2DFiniteElement TriangleFE;
676  const GaussQuad2DFiniteElement TriangleFE;
677  const GaussBiQuad2DFiniteElement QuadrilateralFE;
678 
679 public:
681 
682  virtual const FiniteElement *
684 
685  virtual int DofForGeometry(Geometry::Type GeomType) const;
686 
687  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
688  int Or) const;
689 
690  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
691 };
692 
693 /** Piecewise-cubic discontinuous finite elements in 2D. This class is kept
694  only for backward compatibility, consider using L2_FECollection instead. */
696 {
697 private:
698  const Cubic2DFiniteElement TriangleFE;
699  const BiCubic2DFiniteElement QuadrilateralFE;
700 
701 public:
703 
704  virtual const FiniteElement *
706 
707  virtual int DofForGeometry(Geometry::Type GeomType) const;
708 
709  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
710  int Or) const;
711 
712  virtual const char * Name() const { return "CubicDiscont2D"; }
713 };
714 
715 /** Piecewise-constant discontinuous finite elements in 3D. This class is kept
716  only for backward compatibility, consider using L2_FECollection instead. */
718 {
719 private:
720  const P0TetFiniteElement TetrahedronFE;
721  const P0HexFiniteElement ParallelepipedFE;
722  const L2_WedgeElement WedgeFE;
723 
724 public:
725  Const3DFECollection() : WedgeFE(0) { }
726 
727  virtual const FiniteElement *
729 
730  virtual int DofForGeometry(Geometry::Type GeomType) const;
731 
732  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
733  int Or) const;
734 
735  virtual const char * Name() const { return "Const3D"; }
736 };
737 
738 /** Piecewise-linear discontinuous finite elements in 3D. This class is kept
739  only for backward compatibility, consider using L2_FECollection instead. */
741 {
742 private:
743  const Linear3DFiniteElement TetrahedronFE;
744  const TriLinear3DFiniteElement ParallelepipedFE;
745 
746 public:
748 
749  virtual const FiniteElement *
751 
752  virtual int DofForGeometry(Geometry::Type GeomType) const;
753 
754  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
755  int Or) const;
756 
757  virtual const char * Name() const { return "LinearDiscont3D"; }
758 };
759 
760 /** Piecewise-quadratic discontinuous finite elements in 3D. This class is kept
761  only for backward compatibility, consider using L2_FECollection instead. */
763 {
764 private:
765  const Quadratic3DFiniteElement TetrahedronFE;
766  const LagrangeHexFiniteElement ParallelepipedFE;
767 
768 public:
769  QuadraticDiscont3DFECollection () : ParallelepipedFE(2) { }
770 
771  virtual const FiniteElement *
773 
774  virtual int DofForGeometry(Geometry::Type GeomType) const;
775 
776  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
777  int Or) const;
778 
779  virtual const char * Name() const { return "QuadraticDiscont3D"; }
780 };
781 
782 /// Finite element collection on a macro-element.
784 {
785 private:
786  const PointFiniteElement PointFE;
787  const RefinedLinear1DFiniteElement SegmentFE;
788  const RefinedLinear2DFiniteElement TriangleFE;
789  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
790  const RefinedLinear3DFiniteElement TetrahedronFE;
791  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
792 
793 public:
795 
796  virtual const FiniteElement *
798 
799  virtual int DofForGeometry(Geometry::Type GeomType) const;
800 
801  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
802  int Or) const;
803 
804  virtual const char * Name() const { return "RefinedLinear"; }
805 };
806 
807 /** Lowest order Nedelec finite elements in 3D. This class is kept only for
808  backward compatibility, consider using the new ND_FECollection instead. */
810 {
811 private:
812  const Nedelec1HexFiniteElement HexahedronFE;
813  const Nedelec1TetFiniteElement TetrahedronFE;
814 
815 public:
817 
818  virtual const FiniteElement *
820 
821  virtual int DofForGeometry(Geometry::Type GeomType) const;
822 
823  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
824  int Or) const;
825 
826  virtual const char * Name() const { return "ND1_3D"; }
827 };
828 
829 /** First order Raviart-Thomas finite elements in 3D. This class is kept only
830  for backward compatibility, consider using RT_FECollection instead. */
832 {
833 private:
834  const P0TriangleFiniteElement TriangleFE;
835  const P0QuadFiniteElement QuadrilateralFE;
836  const RT0HexFiniteElement HexahedronFE;
837  const RT0TetFiniteElement TetrahedronFE;
838 public:
840 
841  virtual const FiniteElement *
843 
844  virtual int DofForGeometry(Geometry::Type GeomType) const;
845 
846  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
847  int Or) const;
848 
849  virtual const char * Name() const { return "RT0_3D"; }
850 };
851 
852 /** Second order Raviart-Thomas finite elements in 3D. This class is kept only
853  for backward compatibility, consider using RT_FECollection instead. */
855 {
856 private:
857  const Linear2DFiniteElement TriangleFE;
858  const BiLinear2DFiniteElement QuadrilateralFE;
859  const RT1HexFiniteElement HexahedronFE;
860 public:
862 
863  virtual const FiniteElement *
865 
866  virtual int DofForGeometry(Geometry::Type GeomType) const;
867 
868  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
869  int Or) const;
870 
871  virtual const char * Name() const { return "RT1_3D"; }
872 };
873 
874 /// Discontinuous collection defined locally by a given finite element.
876 {
877 private:
878  char d_name[32];
879  Geometry::Type GeomType;
880  FiniteElement *Local_Element;
881 
882 public:
883  Local_FECollection(const char *fe_name);
884 
886  Geometry::Type _GeomType) const
887  { return (GeomType == _GeomType) ? Local_Element : NULL; }
888  virtual int DofForGeometry(Geometry::Type _GeomType) const
889  { return (GeomType == _GeomType) ? Local_Element->GetDof() : 0; }
890  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
891  int Or) const
892  { return NULL; }
893  virtual const char *Name() const { return d_name; }
894 
895  virtual ~Local_FECollection() { delete Local_Element; }
896 };
897 
898 }
899 
900 #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:838
Abstract class for Finite Elements.
Definition: fe.hpp:232
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:320
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:955
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:289
virtual const char * Name() const
Definition: fe_coll.hpp:438
virtual const char * Name() const
Definition: fe_coll.hpp:779
virtual const char * Name() const
Definition: fe_coll.hpp:690
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1789
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1169
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:1161
For scalar fields; preserves point values.
Definition: fe.hpp:276
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:672
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:197
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1031
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:977
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:255
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:665
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1210
Linear 1D element with nodes 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:1372
Quadratic 1D element with nodes the Gaussian points in [0,1] (trace of RT2)
Definition: fe.hpp:1382
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:1126
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:1240
virtual const char * Name() const
Definition: fe_coll.hpp:804
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:735
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:97
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:2581
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:1975
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:2066
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:574
virtual const char * Name() const
Definition: fe_coll.hpp:592
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:2197
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2615
static const int NumGeom
Definition: geom.hpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:390
const Geometry::Type geom
Definition: ex1.cpp:40
virtual ~L2_FECollection()
Definition: fe_coll.cpp:1991
virtual const char * Name() const
Definition: fe_coll.hpp:629
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2589
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:759
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:616
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:343
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:1719
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:928
void Reset() const
Definition: fe.hpp:2924
Class for quadratic FE on triangle.
Definition: fe.hpp:955
Class for quadratic FE on interval.
Definition: fe.hpp:926
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1377
virtual const char * Name() const
Definition: fe_coll.hpp:222
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2465
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:993
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:662
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:2079
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2215
virtual const char * Name() const
Definition: fe_coll.hpp:334
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:2263
int HasFaceDofs(Geometry::Type GeomType) const
Definition: fe_coll.cpp:25
virtual const char * Name() const
Definition: fe_coll.hpp:505
virtual const char * Name() const
Definition: fe_coll.hpp:712
Finite element collection on a macro-element.
Definition: fe_coll.hpp:783
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:920
Class for refined bi-linear FE on quadrilateral.
Definition: fe.hpp:1511
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:2558
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1184
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:271
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1248
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:626
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:652
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:696
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2231
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1809
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:159
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:415
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:824
virtual const char * Name() const
Definition: fe_coll.hpp:527
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1756
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:879
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:540
Class for bilinear FE on quad with nodes at the 4 Gaussian points.
Definition: fe.hpp:901
Class for refined linear FE on interval.
Definition: fe.hpp:1458
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:90
Class for refined linear FE on triangle.
Definition: fe.hpp:1478
Class for refined linear FE on tetrahedron.
Definition: fe.hpp:1498
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:323
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:196
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:906
Class for constant FE on triangle.
Definition: fe.hpp:1094
virtual int DofForGeometry(Geometry::Type _GeomType) const
Definition: fe_coll.hpp:888
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:293
int GetBasisType() const
Definition: fe_coll.hpp:182
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1365
For scalar fields; preserves volume integrals.
Definition: fe.hpp:277
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:1355
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:1506
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:2497
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2283
Class for refined trilinear FE on a hexahedron.
Definition: fe.hpp:1531
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:875
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:655
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2485
Class for quadratic FE on tetrahedron.
Definition: fe.hpp:1149
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:811
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:963
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:782
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:368
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:1059
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1261
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:394
virtual ~Local_FECollection()
Definition: fe_coll.hpp:895
virtual const char * Name() const
Definition: fe_coll.hpp:614
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:582
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1777
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:1315
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:1277
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:890
Class for linear FE on tetrahedron.
Definition: fe.hpp:1124
Class for linear FE on interval.
Definition: fe.hpp:824
Crouzeix-Raviart finite element on quadrilateral.
Definition: fe.hpp:1198
Class for linear FE on triangle.
Definition: fe.hpp:844
Class for quadratic FE on triangle with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:979
virtual const char * Name() const
Definition: fe_coll.hpp:104
virtual const char * Name() const
Definition: fe_coll.hpp:757
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1737
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:218
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:89
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:463
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:442
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:769
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1298
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1110
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:191
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1001
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:188
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:164
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:303
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:176
H1Ser_FECollection(const int p, const int dim=2)
Definition: fe_coll.hpp:128
virtual const char * Name() const
Definition: fe_coll.hpp:459
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:1438
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:364
virtual const char * Name() const
Definition: fe_coll.hpp:549
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:2568
virtual const char * Name() const
Definition: fe_coll.hpp:482
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:865
virtual const char * Name() const
Definition: fe_coll.hpp:826
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:314
Class for bi-quadratic FE on quadrilateral.
Definition: fe.hpp:995
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:2442
Class for tri-linear FE on cube.
Definition: fe.hpp:1162
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type _GeomType) const
Definition: fe_coll.hpp:885
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:745
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:599
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1285
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1408
virtual const char * Name() const
Definition: fe_coll.hpp:651
Class for bilinear FE on quadrilateral.
Definition: fe.hpp:866
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1339
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:596
virtual const char * Name() const
Definition: fe_coll.hpp:174
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:890
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1322
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:1023
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:731
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1010
virtual const char * Name() const
Definition: fe_coll.hpp:53
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1469
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:267
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:1484
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
Definition: fe.hpp:1034
virtual const char * Name() const
Definition: fe_coll.hpp:570
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:2608
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:639
int dim
Definition: ex24.cpp:43
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1067
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1079
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:215
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1223
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:796
virtual const char * Name() const
Definition: fe_coll.hpp:849
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:679
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1422
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:2243
Tensor products of 1D FEs (only degree 2 is functional)
Definition: fe.hpp:1437
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:557
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:368
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:1200
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:852
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:251
void Reset() const
Definition: fe_coll.hpp:311
Class for cubic FE on tetrahedron.
Definition: fe.hpp:1081
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1044
virtual const char * Name() const
Definition: fe_coll.hpp:893
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
virtual const char * Name() const
Definition: fe_coll.hpp:871
virtual const char * Name() const
Definition: fe_coll.hpp:668
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1095
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2602
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:941
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:894
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:2516
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1147
virtual const char * Name() const
Definition: fe_coll.hpp:411
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:618
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:1393
virtual ~FiniteElementCollection()
Definition: fe_coll.hpp:65
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1134
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:256
Crouzeix-Raviart finite element on triangle.
Definition: fe.hpp:1186
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1456
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:143
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:264