MFEM  v3.4
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, int &fg, int &fo,
40  const int face_info);
41 
42 public:
43  virtual const FiniteElement *
44  FiniteElementForGeometry(int GeomType) const = 0;
45 
46  virtual int DofForGeometry(int GeomType) const = 0;
47 
48  virtual int * DofOrderForOrientation(int GeomType, int Or) const = 0;
49 
50  virtual const char * Name() const { return "Undefined"; }
51 
52  int HasFaceDofs(int GeomType) const;
53 
55  int GeomType) const
56  {
57  return FiniteElementForGeometry(GeomType);
58  }
59 
61 
63 
64  /** @brief Factory method: return a newly allocated FiniteElementCollection
65  according to the given name. */
66  static FiniteElementCollection *New(const char *name);
67 
68  /** @brief Get the local dofs for a given sub-manifold.
69 
70  Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex,
71  1D - edge, 2D - face) including those on its boundary. The local index of
72  the sub-manifold (inside Geom) and its orientation are given by the
73  parameter Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed
74  that 0 <= SDim <= Dim(Geom). */
75  void SubDofOrder(int Geom, int SDim, int Info, Array<int> &dofs) const;
76 };
77 
78 /// Arbitrary order H1-conforming (continuous) finite elements.
80 {
81 
82 protected:
83  int b_type;
84  char h1_name[32];
87  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
88 
89 public:
90  explicit H1_FECollection(const int p, const int dim = 3,
91  const int btype = BasisType::GaussLobatto);
92 
93  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
94  { return H1_Elements[GeomType]; }
95  virtual int DofForGeometry(int GeomType) const
96  { return H1_dof[GeomType]; }
97  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
98  virtual const char *Name() const { return h1_name; }
100 
101  int GetBasisType() const { return b_type; }
102  /// Get the Cartesian to local H1 dof map
103  const int *GetDofMap(int GeomType) const;
104 
105  virtual ~H1_FECollection();
106 };
107 
108 /** Arbitrary order H1-conforming (continuous) finite elements with positive
109  basis functions. */
111 {
112 public:
113  explicit H1Pos_FECollection(const int p, const int dim = 3)
114  : H1_FECollection(p, dim, BasisType::Positive) { }
115 };
116 
117 /** Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the
118  interface between mesh elements (faces,edges,vertices); these are the trace
119  FEs of the H1-conforming FEs. */
121 {
122 public:
123  H1_Trace_FECollection(const int p, const int dim,
124  const int btype = BasisType::GaussLobatto);
125 };
126 
127 /// Arbitrary order "L2-conforming" discontinuous finite elements.
129 {
130 private:
131  int b_type; // BasisType
132  char d_name[32];
135  int *SegDofOrd[2]; // for rotating segment dofs in 1D
136  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
137  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
138 
139 public:
140  L2_FECollection(const int p, const int dim,
141  const int btype = BasisType::GaussLegendre,
142  const int map_type = FiniteElement::VALUE);
143 
144  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
145  { return L2_Elements[GeomType]; }
146  virtual int DofForGeometry(int GeomType) const
147  {
148  if (L2_Elements[GeomType])
149  {
150  return L2_Elements[GeomType]->GetDof();
151  }
152  return 0;
153  }
154  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
155  virtual const char *Name() const { return d_name; }
156 
158  int GeomType) const
159  {
160  return Tr_Elements[GeomType];
161  }
162 
163  int GetBasisType() const { return b_type; }
164 
165  virtual ~L2_FECollection();
166 };
167 
168 /// Declare an alternative name for L2_FECollection = DG_FECollection
170 
171 /// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
173 {
174 protected:
175  int ob_type; // open BasisType
176  char rt_name[32];
179  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
180 
181  // Initialize only the face elements
182  void InitFaces(const int p, const int dim, const int map_type,
183  const bool signs);
184 
185  // Constructor used by the constructor of the RT_Trace_FECollection and
186  // DG_Interface_FECollection classes
187  RT_FECollection(const int p, const int dim, const int map_type,
188  const bool signs,
189  const int ob_type = BasisType::GaussLegendre);
190 
191 public:
192  RT_FECollection(const int p, const int dim,
193  const int cb_type = BasisType::GaussLobatto,
194  const int ob_type = BasisType::GaussLegendre);
195 
196  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
197  { return RT_Elements[GeomType]; }
198  virtual int DofForGeometry(int GeomType) const
199  { return RT_dof[GeomType]; }
200  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
201  virtual const char *Name() const { return rt_name; }
203 
204  virtual ~RT_FECollection();
205 };
206 
207 /** Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the
208  interface between mesh elements (faces); these are the normal trace FEs of
209  the H(div)-conforming FEs. */
211 {
212 public:
213  RT_Trace_FECollection(const int p, const int dim,
214  const int map_type = FiniteElement::INTEGRAL,
215  const int ob_type = BasisType::GaussLegendre);
216 };
217 
218 /** Arbitrary order discontinuous finite elements defined on the interface
219  between mesh elements (faces). The functions in this space are single-valued
220  on each face and are discontinuous across its boundary. */
222 {
223 public:
224  DG_Interface_FECollection(const int p, const int dim,
225  const int map_type = FiniteElement::VALUE,
226  const int ob_type = BasisType::GaussLegendre);
227 };
228 
229 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
231 {
232 protected:
233  char nd_name[32];
236  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
237 
238 public:
239  ND_FECollection(const int p, const int dim,
240  const int cb_type = BasisType::GaussLobatto,
241  const int ob_type = BasisType::GaussLegendre);
242 
243  virtual const FiniteElement *FiniteElementForGeometry(int GeomType) const
244  { return ND_Elements[GeomType]; }
245  virtual int DofForGeometry(int GeomType) const
246  { return ND_dof[GeomType]; }
247  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
248  virtual const char *Name() const { return nd_name; }
250 
251  virtual ~ND_FECollection();
252 };
253 
254 /** Arbitrary order H(curl)-trace finite elements defined on the interface
255  between mesh elements (faces,edges); these are the tangential trace FEs of
256  the H(curl)-conforming FEs. */
258 {
259 public:
260  ND_Trace_FECollection(const int p, const int dim,
261  const int cb_type = BasisType::GaussLobatto,
262  const int ob_type = BasisType::GaussLegendre);
263 };
264 
265 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
267 {
268 private:
269  NURBS1DFiniteElement *SegmentFE;
270  NURBS2DFiniteElement *QuadrilateralFE;
271  NURBS3DFiniteElement *ParallelepipedFE;
272 
273  mutable int mOrder; // >= 1 or VariableOrder
274  // The 'name' can be:
275  // 1) name = "NURBS" + "number", for fixed order, or
276  // 2) name = "NURBS", for VariableOrder.
277  // The name is updated before writing it to a stream, for example, see
278  // FiniteElementSpace::Save().
279  mutable char name[16];
280 
281 public:
282  enum { VariableOrder = -1 };
283 
284  /** @brief The parameter @a Order must be either a positive number, for fixed
285  order, or VariableOrder (default). */
286  explicit NURBSFECollection(int Order = VariableOrder);
287 
288  void Reset() const
289  {
290  SegmentFE->Reset();
291  QuadrilateralFE->Reset();
292  ParallelepipedFE->Reset();
293  }
294 
295  /** @brief Get the order of the NURBS collection: either a positive number,
296  when using fixed order, or VariableOrder. */
297  int GetOrder() const { return mOrder; }
298 
299  /** @brief Set the order and the name, based on the given @a Order: either a
300  positive number for fixed order, or VariableOrder. */
301  void SetOrder(int Order) const;
302 
303  virtual const FiniteElement *
304  FiniteElementForGeometry(int GeomType) const;
305 
306  virtual int DofForGeometry(int GeomType) const;
307 
308  virtual int *DofOrderForOrientation(int GeomType, int Or) const;
309 
310  virtual const char *Name() const { return name; }
311 
313 
314  virtual ~NURBSFECollection();
315 };
316 
317 
318 /// Piecewise-(bi)linear continuous finite elements.
320 {
321 private:
322  const PointFiniteElement PointFE;
323  const Linear1DFiniteElement SegmentFE;
324  const Linear2DFiniteElement TriangleFE;
325  const BiLinear2DFiniteElement QuadrilateralFE;
326  const Linear3DFiniteElement TetrahedronFE;
327  const TriLinear3DFiniteElement ParallelepipedFE;
328 public:
330 
331  virtual const FiniteElement *
332  FiniteElementForGeometry(int GeomType) const;
333 
334  virtual int DofForGeometry(int GeomType) const;
335 
336  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
337 
338  virtual const char * Name() const { return "Linear"; }
339 };
340 
341 /// Piecewise-(bi)quadratic continuous finite elements.
343 {
344 private:
345  const PointFiniteElement PointFE;
346  const Quad1DFiniteElement SegmentFE;
347  const Quad2DFiniteElement TriangleFE;
348  const BiQuad2DFiniteElement QuadrilateralFE;
349  const Quadratic3DFiniteElement TetrahedronFE;
350  const LagrangeHexFiniteElement ParallelepipedFE;
351 
352 public:
353  QuadraticFECollection() : ParallelepipedFE(2) { }
354 
355  virtual const FiniteElement *
356  FiniteElementForGeometry(int GeomType) const;
357 
358  virtual int DofForGeometry(int GeomType) const;
359 
360  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
361 
362  virtual const char * Name() const { return "Quadratic"; }
363 };
364 
365 /// Version of QuadraticFECollection with positive basis functions.
367 {
368 private:
369  const QuadPos1DFiniteElement SegmentFE;
370  const BiQuadPos2DFiniteElement QuadrilateralFE;
371 
372 public:
374 
375  virtual const FiniteElement *
376  FiniteElementForGeometry(int GeomType) const;
377 
378  virtual int DofForGeometry(int GeomType) const;
379 
380  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
381 
382  virtual const char * Name() const { return "QuadraticPos"; }
383 };
384 
385 /// Piecewise-(bi)cubic continuous finite elements.
387 {
388 private:
389  const PointFiniteElement PointFE;
390  const Cubic1DFiniteElement SegmentFE;
391  const Cubic2DFiniteElement TriangleFE;
392  const BiCubic2DFiniteElement QuadrilateralFE;
393  const Cubic3DFiniteElement TetrahedronFE;
394  const LagrangeHexFiniteElement ParallelepipedFE;
395 
396 public:
397  CubicFECollection() : ParallelepipedFE(3) { }
398 
399  virtual const FiniteElement *
400  FiniteElementForGeometry(int GeomType) const;
401 
402  virtual int DofForGeometry(int GeomType) const;
403 
404  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
405 
406  virtual const char * Name() const { return "Cubic"; }
407 };
408 
409 /// Crouzeix-Raviart nonconforming elements in 2D.
411 {
412 private:
413  const P0SegmentFiniteElement SegmentFE;
414  const CrouzeixRaviartFiniteElement TriangleFE;
415  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
416 public:
417  CrouzeixRaviartFECollection() : SegmentFE(1) { }
418 
419  virtual const FiniteElement *
420  FiniteElementForGeometry(int GeomType) const;
421 
422  virtual int DofForGeometry(int GeomType) const;
423 
424  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
425 
426  virtual const char * Name() const { return "CrouzeixRaviart"; }
427 };
428 
429 /// Piecewise-linear nonconforming finite elements in 3D.
431 {
432 private:
433  const P0TriangleFiniteElement TriangleFE;
434  const P1TetNonConfFiniteElement TetrahedronFE;
435  const P0QuadFiniteElement QuadrilateralFE;
436  const RotTriLinearHexFiniteElement ParallelepipedFE;
437 
438 public:
440 
441  virtual const FiniteElement *
442  FiniteElementForGeometry(int GeomType) const;
443 
444  virtual int DofForGeometry(int GeomType) const;
445 
446  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
447 
448  virtual const char * Name() const { return "LinearNonConf3D"; }
449 };
450 
451 
452 /** First order Raviart-Thomas finite elements in 2D. This class is kept only
453  for backward compatibility, consider using RT_FECollection instead. */
455 {
456 private:
457  const P0SegmentFiniteElement SegmentFE; // normal component on edge
458  const RT0TriangleFiniteElement TriangleFE;
459  const RT0QuadFiniteElement QuadrilateralFE;
460 public:
461  RT0_2DFECollection() : SegmentFE(0) { }
462 
463  virtual const FiniteElement *
464  FiniteElementForGeometry(int GeomType) const;
465 
466  virtual int DofForGeometry(int GeomType) const;
467 
468  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
469 
470  virtual const char * Name() const { return "RT0_2D"; }
471 };
472 
473 /** Second order Raviart-Thomas finite elements in 2D. This class is kept only
474  for backward compatibility, consider using RT_FECollection instead. */
476 {
477 private:
478  const P1SegmentFiniteElement SegmentFE; // normal component on edge
479  const RT1TriangleFiniteElement TriangleFE;
480  const RT1QuadFiniteElement QuadrilateralFE;
481 public:
483 
484  virtual const FiniteElement *
485  FiniteElementForGeometry(int GeomType) const;
486 
487  virtual int DofForGeometry(int GeomType) const;
488 
489  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
490 
491  virtual const char * Name() const { return "RT1_2D"; }
492 };
493 
494 /** Third order Raviart-Thomas finite elements in 2D. This class is kept only
495  for backward compatibility, consider using RT_FECollection instead. */
497 {
498 private:
499  const P2SegmentFiniteElement SegmentFE; // normal component on edge
500  const RT2TriangleFiniteElement TriangleFE;
501  const RT2QuadFiniteElement QuadrilateralFE;
502 public:
504 
505  virtual const FiniteElement *
506  FiniteElementForGeometry(int GeomType) const;
507 
508  virtual int DofForGeometry(int GeomType) const;
509 
510  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
511 
512  virtual const char * Name() const { return "RT2_2D"; }
513 };
514 
515 /** Piecewise-constant discontinuous finite elements in 2D. This class is kept
516  only for backward compatibility, consider using L2_FECollection instead. */
518 {
519 private:
520  const P0TriangleFiniteElement TriangleFE;
521  const P0QuadFiniteElement QuadrilateralFE;
522 public:
524 
525  virtual const FiniteElement *
526  FiniteElementForGeometry(int GeomType) const;
527 
528  virtual int DofForGeometry(int GeomType) const;
529 
530  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
531 
532  virtual const char * Name() const { return "Const2D"; }
533 };
534 
535 /** Piecewise-linear discontinuous finite elements in 2D. This class is kept
536  only for backward compatibility, consider using L2_FECollection instead. */
538 {
539 private:
540  const Linear2DFiniteElement TriangleFE;
541  const BiLinear2DFiniteElement QuadrilateralFE;
542 
543 public:
545 
546  virtual const FiniteElement *
547  FiniteElementForGeometry(int GeomType) const;
548 
549  virtual int DofForGeometry(int GeomType) const;
550 
551  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
552 
553  virtual const char * Name() const { return "LinearDiscont2D"; }
554 };
555 
556 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
558 {
559 private:
560  // const CrouzeixRaviartFiniteElement TriangleFE;
561  const GaussLinear2DFiniteElement TriangleFE;
562  const GaussBiLinear2DFiniteElement QuadrilateralFE;
563 
564 public:
566 
567  virtual const FiniteElement *
568  FiniteElementForGeometry(int GeomType) const;
569 
570  virtual int DofForGeometry(int GeomType) const;
571 
572  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
573 
574  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
575 };
576 
577 /// Linear (P1) finite elements on quadrilaterals.
579 {
580 private:
581  const P1OnQuadFiniteElement QuadrilateralFE;
582 public:
584  virtual const FiniteElement *
585  FiniteElementForGeometry(int GeomType) const;
586  virtual int DofForGeometry(int GeomType) const;
587  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
588  virtual const char * Name() const { return "P1OnQuad"; }
589 };
590 
591 /** Piecewise-quadratic discontinuous finite elements in 2D. This class is kept
592  only for backward compatibility, consider using L2_FECollection instead. */
594 {
595 private:
596  const Quad2DFiniteElement TriangleFE;
597  const BiQuad2DFiniteElement QuadrilateralFE;
598 
599 public:
601 
602  virtual const FiniteElement *
603  FiniteElementForGeometry(int GeomType) const;
604 
605  virtual int DofForGeometry(int GeomType) const;
606 
607  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
608 
609  virtual const char * Name() const { return "QuadraticDiscont2D"; }
610 };
611 
612 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
614 {
615 private:
616  const BiQuadPos2DFiniteElement QuadrilateralFE;
617 
618 public:
620  virtual const FiniteElement *
621  FiniteElementForGeometry(int GeomType) const;
622  virtual int DofForGeometry(int GeomType) const;
623  virtual int * DofOrderForOrientation(int GeomType, int Or) const
624  { return NULL; }
625  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
626 };
627 
628 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
630 {
631 private:
632  // const Quad2DFiniteElement TriangleFE;
633  const GaussQuad2DFiniteElement TriangleFE;
634  const GaussBiQuad2DFiniteElement QuadrilateralFE;
635 
636 public:
638 
639  virtual const FiniteElement *
640  FiniteElementForGeometry(int GeomType) const;
641 
642  virtual int DofForGeometry(int GeomType) const;
643 
644  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
645 
646  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
647 };
648 
649 /** Piecewise-cubic discontinuous finite elements in 2D. This class is kept
650  only for backward compatibility, consider using L2_FECollection instead. */
652 {
653 private:
654  const Cubic2DFiniteElement TriangleFE;
655  const BiCubic2DFiniteElement QuadrilateralFE;
656 
657 public:
659 
660  virtual const FiniteElement *
661  FiniteElementForGeometry(int GeomType) const;
662 
663  virtual int DofForGeometry(int GeomType) const;
664 
665  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
666 
667  virtual const char * Name() const { return "CubicDiscont2D"; }
668 };
669 
670 /** Piecewise-constant discontinuous finite elements in 3D. This class is kept
671  only for backward compatibility, consider using L2_FECollection instead. */
673 {
674 private:
675  const P0TetFiniteElement TetrahedronFE;
676  const P0HexFiniteElement ParallelepipedFE;
677 
678 public:
680 
681  virtual const FiniteElement *
682  FiniteElementForGeometry(int GeomType) const;
683 
684  virtual int DofForGeometry(int GeomType) const;
685 
686  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
687 
688  virtual const char * Name() const { return "Const3D"; }
689 };
690 
691 /** Piecewise-linear discontinuous finite elements in 3D. This class is kept
692  only for backward compatibility, consider using L2_FECollection instead. */
694 {
695 private:
696  const Linear3DFiniteElement TetrahedronFE;
697  const TriLinear3DFiniteElement ParallelepipedFE;
698 
699 public:
701 
702  virtual const FiniteElement *
703  FiniteElementForGeometry(int GeomType) const;
704 
705  virtual int DofForGeometry(int GeomType) const;
706 
707  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
708 
709  virtual const char * Name() const { return "LinearDiscont3D"; }
710 };
711 
712 /** Piecewise-quadratic discontinuous finite elements in 3D. This class is kept
713  only for backward compatibility, consider using L2_FECollection instead. */
715 {
716 private:
717  const Quadratic3DFiniteElement TetrahedronFE;
718  const LagrangeHexFiniteElement ParallelepipedFE;
719 
720 public:
721  QuadraticDiscont3DFECollection () : ParallelepipedFE(2) { }
722 
723  virtual const FiniteElement *
724  FiniteElementForGeometry(int GeomType) const;
725 
726  virtual int DofForGeometry(int GeomType) const;
727 
728  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
729 
730  virtual const char * Name() const { return "QuadraticDiscont3D"; }
731 };
732 
733 /// Finite element collection on a macro-element.
735 {
736 private:
737  const PointFiniteElement PointFE;
738  const RefinedLinear1DFiniteElement SegmentFE;
739  const RefinedLinear2DFiniteElement TriangleFE;
740  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
741  const RefinedLinear3DFiniteElement TetrahedronFE;
742  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
743 
744 public:
746 
747  virtual const FiniteElement *
748  FiniteElementForGeometry(int GeomType) const;
749 
750  virtual int DofForGeometry(int GeomType) const;
751 
752  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
753 
754  virtual const char * Name() const { return "RefinedLinear"; }
755 };
756 
757 /** Lowest order Nedelec finite elements in 3D. This class is kept only for
758  backward compatibility, consider using the new ND_FECollection instead. */
760 {
761 private:
762  const Nedelec1HexFiniteElement HexahedronFE;
763  const Nedelec1TetFiniteElement TetrahedronFE;
764 
765 public:
767 
768  virtual const FiniteElement *
769  FiniteElementForGeometry(int GeomType) const;
770 
771  virtual int DofForGeometry(int GeomType) const;
772 
773  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
774 
775  virtual const char * Name() const { return "ND1_3D"; }
776 };
777 
778 /** First order Raviart-Thomas finite elements in 3D. This class is kept only
779  for backward compatibility, consider using RT_FECollection instead. */
781 {
782 private:
783  const P0TriangleFiniteElement TriangleFE;
784  const P0QuadFiniteElement QuadrilateralFE;
785  const RT0HexFiniteElement HexahedronFE;
786  const RT0TetFiniteElement TetrahedronFE;
787 public:
789 
790  virtual const FiniteElement *
791  FiniteElementForGeometry(int GeomType) const;
792 
793  virtual int DofForGeometry(int GeomType) const;
794 
795  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
796 
797  virtual const char * Name() const { return "RT0_3D"; }
798 };
799 
800 /** Second order Raviart-Thomas finite elements in 3D. This class is kept only
801  for backward compatibility, consider using RT_FECollection instead. */
803 {
804 private:
805  const Linear2DFiniteElement TriangleFE;
806  const BiLinear2DFiniteElement QuadrilateralFE;
807  const RT1HexFiniteElement HexahedronFE;
808 public:
810 
811  virtual const FiniteElement *
812  FiniteElementForGeometry(int GeomType) const;
813 
814  virtual int DofForGeometry(int GeomType) const;
815 
816  virtual int * DofOrderForOrientation(int GeomType, int Or) const;
817 
818  virtual const char * Name() const { return "RT1_3D"; }
819 };
820 
821 /// Discontinuous collection defined locally by a given finite element.
823 {
824 private:
825  char d_name[32];
826  int GeomType;
827  FiniteElement *Local_Element;
828 
829 public:
830  Local_FECollection(const char *fe_name);
831 
832  virtual const FiniteElement *FiniteElementForGeometry(int _GeomType) const
833  { return (GeomType == _GeomType) ? Local_Element : NULL; }
834  virtual int DofForGeometry(int _GeomType) const
835  { return (GeomType == _GeomType) ? Local_Element->GetDof() : 0; }
836  virtual int *DofOrderForOrientation(int GeomType, int Or) const
837  { return NULL; }
838  virtual const char *Name() const { return d_name; }
839 
840  virtual ~Local_FECollection() { delete Local_Element; }
841 };
842 
843 }
844 
845 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:140
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:297
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:266
virtual const char * Name() const
Definition: fe_coll.hpp:406
virtual const char * Name() const
Definition: fe_coll.hpp:730
virtual const char * Name() const
Definition: fe_coll.hpp:646
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1688
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1251
For scalar fields; preserves point values.
Definition: fe.hpp:180
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:629
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:178
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:234
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1356
Linear 1D element with nodes 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:1208
Quadratic 1D element with nodes the Gaussian points in [0,1] (trace of RT2)
Definition: fe.hpp:1218
virtual const char * Name() const
Definition: fe_coll.hpp:754
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:38
virtual const char * Name() const
Definition: fe_coll.hpp:688
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:2452
void SubDofOrder(int Geom, int SDim, int Info, Array< int > &dofs) const
Get the local dofs for a given sub-manifold.
Definition: fe_coll.cpp:361
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:1941
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:623
virtual const char * Name() const
Definition: fe_coll.hpp:553
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2485
static const int NumGeom
Definition: geom.hpp:34
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2479
virtual const char * Name() const
Definition: fe_coll.hpp:362
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1419
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1179
const Geometry::Type geom
Definition: ex1.cpp:40
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1097
virtual ~L2_FECollection()
Definition: fe_coll.cpp:1868
virtual const char * Name() const
Definition: fe_coll.hpp:588
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:712
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1340
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:942
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1000
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:978
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:908
Piecewise-(bi)linear continuous finite elements.
Definition: fe_coll.hpp:319
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1447
void Reset() const
Definition: fe.hpp:2549
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1279
Class for quadratic FE on triangle.
Definition: fe.hpp:792
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:818
Class for quadratic FE on interval.
Definition: fe.hpp:763
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1286
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1385
virtual const char * Name() const
Definition: fe_coll.hpp:201
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2336
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:198
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1131
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1117
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:900
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:1954
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2089
virtual const char * Name() const
Definition: fe_coll.hpp:310
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1318
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:2137
virtual const char * Name() const
Definition: fe_coll.hpp:470
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1215
virtual const char * Name() const
Definition: fe_coll.hpp:667
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:845
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1068
Finite element collection on a macro-element.
Definition: fe_coll.hpp:734
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:791
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:547
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:570
Class for refined bi-linear FE on quadrilateral.
Definition: fe.hpp:1347
Possible basis types. Note that not all elements can use all BasisType(s).
Definition: fe.hpp:27
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:243
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:832
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default)...
Definition: fe_coll.cpp:2429
virtual const char * Name() const
Definition: fe_coll.hpp:248
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:970
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2105
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1708
virtual int DofForGeometry(int _GeomType) const
Definition: fe_coll.hpp:834
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:886
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2072
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1401
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:386
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:934
virtual const char * Name() const
Definition: fe_coll.hpp:491
int dim
Definition: ex3.cpp:47
Class for bilinear FE on quad with nodes at the 4 Gaussian points.
Definition: fe.hpp:738
Class for refined linear FE on interval.
Definition: fe.hpp:1294
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1432
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:86
Class for refined linear FE on triangle.
Definition: fe.hpp:1314
Class for refined linear FE on tetrahedron.
Definition: fe.hpp:1334
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:177
Class for constant FE on triangle.
Definition: fe.hpp:931
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:286
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:636
int GetBasisType() const
Definition: fe_coll.hpp:163
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:920
For scalar fields; preserves volume integrals.
Definition: fe.hpp:181
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:563
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:987
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1191
int GetBasisType() const
Definition: fe_coll.hpp:101
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1469
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:2368
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2157
Class for refined trilinear FE on a hexahedron.
Definition: fe.hpp:1367
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:822
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:611
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:613
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2356
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1207
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:749
Class for quadratic FE on tetrahedron.
Definition: fe.hpp:986
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:836
const int * GetDofMap(int GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1655
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1020
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:874
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:366
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1139
virtual ~Local_FECollection()
Definition: fe_coll.hpp:840
virtual const char * Name() const
Definition: fe_coll.hpp:574
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1302
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1676
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:531
Class for linear FE on tetrahedron.
Definition: fe.hpp:961
Class for linear FE on interval.
Definition: fe.hpp:661
Crouzeix-Raviart finite element on quadrilateral.
Definition: fe.hpp:1034
Class for linear FE on triangle.
Definition: fe.hpp:681
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1371
Class for quadratic FE on triangle with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:816
virtual const char * Name() const
Definition: fe_coll.hpp:98
virtual const char * Name() const
Definition: fe_coll.hpp:709
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1636
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:85
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:430
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:678
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:410
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:172
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:169
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:739
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:113
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:725
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:296
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:804
virtual const char * Name() const
Definition: fe_coll.hpp:426
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:146
virtual const char * Name() const
Definition: fe_coll.hpp:338
virtual const char * Name() const
Definition: fe_coll.hpp:512
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:157
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1227
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:2439
virtual const char * Name() const
Definition: fe_coll.hpp:448
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1034
virtual const char * Name() const
Definition: fe_coll.hpp:775
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:217
Class for bi-quadratic FE on quadrilateral.
Definition: fe.hpp:832
Class for tri-linear FE on cube.
Definition: fe.hpp:999
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1042
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:602
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:955
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:196
virtual const char * Name() const
Definition: fe_coll.hpp:609
Class for bilinear FE on quadrilateral.
Definition: fe.hpp:703
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1263
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1105
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1082
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:557
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:2460
virtual int * DofOrderForOrientation(int GeomType, int Or) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:155
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:44
virtual const FiniteElement * TraceFiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:54
Class for linear FE on triangle with nodes at the 3 &quot;Gaussian&quot; points.
Definition: fe.hpp:727
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1243
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1328
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:245
virtual const char * Name() const
Definition: fe_coll.hpp:50
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:586
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:93
Bi-quadratic element on quad with nodes at the 9 Gaussian points.
Definition: fe.hpp:871
virtual const char * Name() const
Definition: fe_coll.hpp:532
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:859
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:2313
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1619
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:776
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.hpp:95
virtual const char * Name() const
Definition: fe_coll.hpp:797
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:2117
Tensor products of 1D FEs (only degree 2 is functional)
Definition: fe.hpp:1273
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:342
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.hpp:144
int HasFaceDofs(int GeomType) const
Definition: fe_coll.cpp:25
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:230
void Reset() const
Definition: fe_coll.hpp:288
Class for cubic FE on tetrahedron.
Definition: fe.hpp:918
static void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo, int &nf, int &f, int &fg, int &fo, const int face_info)
Definition: fe_coll.cpp:316
virtual const char * Name() const
Definition: fe_coll.hpp:838
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1169
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.cpp:1853
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
virtual const char * Name() const
Definition: fe_coll.hpp:818
virtual const char * Name() const
Definition: fe_coll.hpp:625
virtual const FiniteElement * FiniteElementForGeometry(int _GeomType) const
Definition: fe_coll.hpp:832
virtual int DofForGeometry(int GeomType) const =0
virtual int * DofOrderForOrientation(int GeomType, int Or) const
Definition: fe_coll.hpp:623
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:2473
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:2387
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1053
virtual const char * Name() const
Definition: fe_coll.hpp:382
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:1153
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:578
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:662
virtual ~FiniteElementCollection()
Definition: fe_coll.hpp:62
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:235
virtual int DofForGeometry(int GeomType) const
Definition: fe_coll.cpp:762
Crouzeix-Raviart finite element on triangle.
Definition: fe.hpp:1022
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:646
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const =0
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:128
virtual const FiniteElement * FiniteElementForGeometry(int GeomType) const
Definition: fe_coll.cpp:1008