MFEM  v4.2.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 /** @brief Collection of finite elements from the same family in multiple
23  dimensions. This class is used to match the degrees of freedom of a
24  FiniteElementSpace between elements, and to provide the finite element
25  restriction from an 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  /** @brief Enumeration for ContType: defines the continuity of the field
44  across element interfaces. */
45  enum { CONTINUOUS, ///< Field is continuous across element interfaces
46  TANGENTIAL, ///< Tangential components of vector field
47  NORMAL, ///< Normal component of vector field
48  DISCONTINUOUS ///< Field is discontinuous across element interfaces
49  };
50 
51  virtual const FiniteElement *
52  FiniteElementForGeometry(Geometry::Type GeomType) const = 0;
53 
54  virtual int DofForGeometry(Geometry::Type GeomType) const = 0;
55 
56  /** @brief Returns an array, say p, that maps a local permuted index i to
57  a local base index: base_i = p[i]. */
58  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
59  int Or) const = 0;
60 
61  virtual const char * Name() const { return "Undefined"; }
62 
63  virtual int GetContType() const = 0;
64 
65  int HasFaceDofs(Geometry::Type GeomType) const;
66 
68  Geometry::Type GeomType) const
69  {
70  return FiniteElementForGeometry(GeomType);
71  }
72 
74 
76 
77  /** @brief Factory method: return a newly allocated FiniteElementCollection
78  according to the given name. */
79  /**
80  | FEC Name | Space | Order | BasisType | FiniteElement::MapT | Notes |
81  | :------: | :---: | :---: | :-------: | :-----: | :---: |
82  | H1_[DIM]_[ORDER] | H1 | * | 1 | VALUE | H1 nodal elements |
83  | H1@[BTYPE]_[DIM]_[ORDER] | H1 | * | * | VALUE | H1 nodal elements |
84  | H1Pos_[DIM]_[ORDER] | H1 | * | 1 | VALUE | H1 nodal elements |
85  | H1Pos_Trace_[DIM]_[ORDER] | H^{1/2} | * | 2 | VALUE | H^{1/2}-conforming trace elements for H1 defined on the interface between mesh elements (faces,edges,vertices) |
86  | H1_Trace_[DIM]_[ORDER] | H^{1/2} | * | 1 | VALUE | H^{1/2}-conforming trace elements for H1 defined on the interface between mesh elements (faces,edges,vertices) |
87  | H1_Trace@[BTYPE]_[DIM]_[ORDER] | H^{1/2} | * | 1 | VALUE | H^{1/2}-conforming trace elements for H1 defined on the interface between mesh elements (faces,edges,vertices) |
88  | ND_[DIM]_[ORDER] | H(curl) | * | 1 / 0 | H_CURL | Nedelec vector elements |
89  | ND@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(curl) | * | * / * | H_CURL | Nedelec vector elements |
90  | ND_Trace_[DIM]_[ORDER] | H^{1/2} | * | 1 / 0 | H_CURL | H^{1/2}-conforming trace elements for H(curl) defined on the interface between mesh elements (faces) |
91  | ND_Trace@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H^{1/2} | * | 1 / 0 | H_CURL | H^{1/2}-conforming trace elements for H(curl) defined on the interface between mesh elements (faces) |
92  | RT_[DIM]_[ORDER] | H(div) | * | 1 / 0 | H_DIV | Raviart-Thomas vector elements |
93  | RT@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(div) | * | * / * | H_DIV | Raviart-Thomas vector elements |
94  | RT_Trace_[DIM]_[ORDER] | H^{1/2} | * | 1 / 0 | INTEGRAL | H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces) |
95  | RT_ValTrace_[DIM]_[ORDER] | H^{1/2} | * | 1 / 0 | VALUE | H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces) |
96  | RT_Trace@[BTYPE]_[DIM]_[ORDER] | H^{1/2} | * | 1 / 0 | INTEGRAL | H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces) |
97  | RT_ValTrace@[BTYPE]_[DIM]_[ORDER] | H^{1/2} | * | 1 / 0 | VALUE | H^{1/2}-conforming trace elements for H(div) defined on the interface between mesh elements (faces) |
98  | L2_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinous L2 elements |
99  | L2_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinous L2 elements |
100  | L2Int_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinous L2 elements |
101  | L2Int_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinous L2 elements |
102  | DG_Iface_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
103  | DG_Iface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
104  | DG_IntIface_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
105  | DG_IntIface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
106  | NURBS[ORDER] | - | * | - | VALUE | Non-Uniform Rational B-Splines (NURBS) elements |
107  | LinearNonConf3D | - | 1 | 1 | VALUE | Piecewise-linear nonconforming finite elements in 3D |
108  | CrouzeixRaviart | - | - | - | - | Crouzeix-Raviart nonconforming elements in 2D |
109  | Local_[FENAME] | - | - | - | - | Special collection that builds a local version out of the FENAME collection |
110  |-|-|-|-|-|-|
111  | Linear | H1 | 1 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
112  | Quadratic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
113  | QuadraticPos | H1 | 2 | 2 | VALUE | Left in for backward compatibility, consider using H1_ |
114  | Cubic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
115  | Const2D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
116  | Const3D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
117  | LinearDiscont2D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
118  | GaussLinearDiscont2D | L2 | 1 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
119  | P1OnQuad | H1 | 1 | 1 | VALUE | Linear P1 element with 3 nodes on a square |
120  | QuadraticDiscont2D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
121  | QuadraticPosDiscont2D | L2 | 2 | 2 | VALUE | Left in for backward compatibility, consider using L2_ |
122  | GaussQuadraticDiscont2D | L2 | 2 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
123  | CubicDiscont2D | L2 | 3 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
124  | LinearDiscont3D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
125  | QuadraticDiscont3D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
126  | ND1_3D | H(Curl) | 1 | 1 / 0 | H_CURL | Left in for backward compatibility, consider using ND_ |
127  | RT0_2D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
128  | RT1_2D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
129  | RT2_2D | H(Div) | 3 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
130  | RT0_3D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
131  | RT1_3D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
132 
133  | Tag | Description |
134  | :------: | :--------: |
135  | [DIM] | Dimension of the elements (1D, 2D, 3D) |
136  | [ORDER] | Approximation order of the elements (P0, P1, P2, ...) |
137  | [BTYPE] | BasisType of the element (0-GaussLegendre, 1 - GaussLobatto, 2-Bernstein, 3-OpenUniform, 4-CloseUniform, 5-OpenHalfUniform) |
138  | [OBTYPE] | Open BasisType of the element for elements which have both types |
139  | [CBTYPE] | Closed BasisType of the element for elements which have both types |
140 
141  [FENAME] Is a special case for the Local FEC which generates a local version of a given
142  FEC. It is selected from one of (BiCubic2DFiniteElement, Quad_Q3, Nedelec1HexFiniteElement,
143  Hex_ND1, H1_[DIM]_[ORDER],H1Pos_[DIM]_[ORDER], L2_[DIM]_[ORDER] )
144  */
145  static FiniteElementCollection *New(const char *name);
146 
147  /** @brief Get the local dofs for a given sub-manifold.
148 
149  Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex, 1D
150  - edge, 2D - face) including those on its boundary. The local index of the
151  sub-manifold (inside Geom) and its orientation are given by the parameter
152  Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed that 0 <=
153  SDim <= Dim(Geom). */
154  void SubDofOrder(Geometry::Type Geom, int SDim, int Info,
155  Array<int> &dofs) const;
156 };
157 
158 /// Arbitrary order H1-conforming (continuous) finite elements.
160 {
161 
162 protected:
163  int b_type;
164  char h1_name[32];
167  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8], *TetDofOrd[24];
168 
169 public:
170  explicit H1_FECollection(const int p, const int dim = 3,
171  const int btype = BasisType::GaussLobatto);
172 
174  Geometry::Type GeomType) const
175  { return H1_Elements[GeomType]; }
176  virtual int DofForGeometry(Geometry::Type GeomType) const
177  { return H1_dof[GeomType]; }
178  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
179  int Or) const;
180  virtual const char *Name() const { return h1_name; }
181  virtual int GetContType() const { return CONTINUOUS; }
183 
184  int GetBasisType() const { return b_type; }
185  /// Get the Cartesian to local H1 dof map
186  const int *GetDofMap(Geometry::Type GeomType) const;
187 
188  virtual ~H1_FECollection();
189 };
190 
191 /** @brief Arbitrary order H1-conforming (continuous) finite elements with
192  positive basis functions. */
194 {
195 public:
196  explicit H1Pos_FECollection(const int p, const int dim = 3)
197  : H1_FECollection(p, dim, BasisType::Positive) { }
198 };
199 
200 
201 /** Arbitrary order H1-conforming (continuous) serendipity finite elements;
202  Current implementation works in 2D only; 3D version is in development. */
204 {
205 public:
206  explicit H1Ser_FECollection(const int p, const int dim = 2)
207  : H1_FECollection(p, dim, BasisType::Serendipity) { };
208 };
209 
210 /** @brief Arbitrary order "H^{1/2}-conforming" trace finite elements defined on
211  the interface between mesh elements (faces,edges,vertices); these are the
212  trace FEs of the H1-conforming FEs. */
214 {
215 public:
216  H1_Trace_FECollection(const int p, const int dim,
217  const int btype = BasisType::GaussLobatto);
218 };
219 
220 /// Arbitrary order "L2-conforming" discontinuous finite elements.
222 {
223 private:
224  int b_type; // BasisType
225  char d_name[32];
228  int *SegDofOrd[2]; // for rotating segment dofs in 1D
229  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
230  int *TetDofOrd[24]; // for rotating tetrahedron dofs in 3D
231  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
232 
233 public:
234  L2_FECollection(const int p, const int dim,
235  const int btype = BasisType::GaussLegendre,
236  const int map_type = FiniteElement::VALUE);
237 
239  Geometry::Type GeomType) const
240  {
241  return L2_Elements[GeomType];
242  }
243  virtual int DofForGeometry(Geometry::Type GeomType) const
244  {
245  if (L2_Elements[GeomType])
246  {
247  return L2_Elements[GeomType]->GetDof();
248  }
249  return 0;
250  }
251  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
252  int Or) const;
253  virtual const char *Name() const { return d_name; }
254 
255  virtual int GetContType() const { return DISCONTINUOUS; }
256 
258  Geometry::Type GeomType) const
259  {
260  return Tr_Elements[GeomType];
261  }
262 
263  int GetBasisType() const { return b_type; }
264 
265  virtual ~L2_FECollection();
266 };
267 
268 /// Declare an alternative name for L2_FECollection = DG_FECollection
270 
271 /// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
273 {
274 protected:
275  int ob_type; // open BasisType
276  char rt_name[32];
279  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
280 
281  // Initialize only the face elements
282  void InitFaces(const int p, const int dim, const int map_type,
283  const bool signs);
284 
285  // Constructor used by the constructor of the RT_Trace_FECollection and
286  // DG_Interface_FECollection classes
287  RT_FECollection(const int p, const int dim, const int map_type,
288  const bool signs,
289  const int ob_type = BasisType::GaussLegendre);
290 
291 public:
292  RT_FECollection(const int p, const int dim,
293  const int cb_type = BasisType::GaussLobatto,
294  const int ob_type = BasisType::GaussLegendre);
295 
297  Geometry::Type GeomType) const
298  { return RT_Elements[GeomType]; }
299  virtual int DofForGeometry(Geometry::Type GeomType) const
300  { return RT_dof[GeomType]; }
301  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
302  int Or) const;
303  virtual const char *Name() const { return rt_name; }
304  virtual int GetContType() const { return NORMAL; }
306 
307  virtual ~RT_FECollection();
308 };
309 
310 /** @brief Arbitrary order "H^{-1/2}-conforming" face finite elements defined on
311  the interface between mesh elements (faces); these are the normal trace FEs
312  of the H(div)-conforming FEs. */
314 {
315 public:
316  RT_Trace_FECollection(const int p, const int dim,
317  const int map_type = FiniteElement::INTEGRAL,
318  const int ob_type = BasisType::GaussLegendre);
319 };
320 
321 /** Arbitrary order discontinuous finite elements defined on the interface
322  between mesh elements (faces). The functions in this space are single-valued
323  on each face and are discontinuous across its boundary. */
325 {
326 public:
327  DG_Interface_FECollection(const int p, const int dim,
328  const int map_type = FiniteElement::VALUE,
329  const int ob_type = BasisType::GaussLegendre);
330 };
331 
332 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
334 {
335 protected:
336  char nd_name[32];
339  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
340 
341 public:
342  ND_FECollection(const int p, const int dim,
343  const int cb_type = BasisType::GaussLobatto,
344  const int ob_type = BasisType::GaussLegendre);
345 
347  const
348  { return ND_Elements[GeomType]; }
349  virtual int DofForGeometry(Geometry::Type GeomType) const
350  { return ND_dof[GeomType]; }
351  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
352  int Or) const;
353  virtual const char *Name() const { return nd_name; }
354  virtual int GetContType() const { return TANGENTIAL; }
356 
357  virtual ~ND_FECollection();
358 };
359 
360 /** @brief Arbitrary order H(curl)-trace finite elements defined on the
361  interface between mesh elements (faces,edges); these are the tangential
362  trace FEs of the H(curl)-conforming FEs. */
364 {
365 public:
366  ND_Trace_FECollection(const int p, const int dim,
367  const int cb_type = BasisType::GaussLobatto,
368  const int ob_type = BasisType::GaussLegendre);
369 };
370 
371 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
373 {
374 private:
375  NURBS1DFiniteElement *SegmentFE;
376  NURBS2DFiniteElement *QuadrilateralFE;
377  NURBS3DFiniteElement *ParallelepipedFE;
378 
379  mutable int mOrder; // >= 1 or VariableOrder
380  // The 'name' can be:
381  // 1) name = "NURBS" + "number", for fixed order, or
382  // 2) name = "NURBS", for VariableOrder.
383  // The name is updated before writing it to a stream, for example, see
384  // FiniteElementSpace::Save().
385  mutable char name[16];
386 
387 public:
388  enum { VariableOrder = -1 };
389 
390  /** @brief The parameter @a Order must be either a positive number, for fixed
391  order, or VariableOrder (default). */
392  explicit NURBSFECollection(int Order = VariableOrder);
393 
394  void Reset() const
395  {
396  SegmentFE->Reset();
397  QuadrilateralFE->Reset();
398  ParallelepipedFE->Reset();
399  }
400 
401  /** @brief Get the order of the NURBS collection: either a positive number,
402  when using fixed order, or VariableOrder. */
403  int GetOrder() const { return mOrder; }
404 
405  /** @brief Set the order and the name, based on the given @a Order: either a
406  positive number for fixed order, or VariableOrder. */
407  void SetOrder(int Order) const;
408 
409  virtual const FiniteElement *
411 
412  virtual int DofForGeometry(Geometry::Type GeomType) const;
413 
414  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
415  int Or) const;
416 
417  virtual const char *Name() const { return name; }
418 
419  virtual int GetContType() const { return CONTINUOUS; }
420 
422 
423  virtual ~NURBSFECollection();
424 };
425 
426 
427 /// Piecewise-(bi/tri)linear continuous finite elements.
429 {
430 private:
431  const PointFiniteElement PointFE;
432  const Linear1DFiniteElement SegmentFE;
433  const Linear2DFiniteElement TriangleFE;
434  const BiLinear2DFiniteElement QuadrilateralFE;
435  const Linear3DFiniteElement TetrahedronFE;
436  const TriLinear3DFiniteElement ParallelepipedFE;
437  const H1_WedgeElement WedgeFE;
438 public:
439  LinearFECollection() : WedgeFE(1) { }
440 
441  virtual const FiniteElement *
443 
444  virtual int DofForGeometry(Geometry::Type GeomType) const;
445 
446  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
447  int Or) const;
448 
449  virtual const char * Name() const { return "Linear"; }
450 
451  virtual int GetContType() const { return CONTINUOUS; }
452 };
453 
454 /// Piecewise-(bi)quadratic continuous finite elements.
456 {
457 private:
458  const PointFiniteElement PointFE;
459  const Quad1DFiniteElement SegmentFE;
460  const Quad2DFiniteElement TriangleFE;
461  const BiQuad2DFiniteElement QuadrilateralFE;
462  const Quadratic3DFiniteElement TetrahedronFE;
463  const LagrangeHexFiniteElement ParallelepipedFE;
464  const H1_WedgeElement WedgeFE;
465 
466 public:
467  QuadraticFECollection() : ParallelepipedFE(2), WedgeFE(2) { }
468 
469  virtual const FiniteElement *
471 
472  virtual int DofForGeometry(Geometry::Type GeomType) const;
473 
474  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
475  int Or) const;
476 
477  virtual const char * Name() const { return "Quadratic"; }
478 
479  virtual int GetContType() const { return CONTINUOUS; }
480 };
481 
482 /// Version of QuadraticFECollection with positive basis functions.
484 {
485 private:
486  const QuadPos1DFiniteElement SegmentFE;
487  const BiQuadPos2DFiniteElement QuadrilateralFE;
488 
489 public:
491 
492  virtual const FiniteElement *
494 
495  virtual int DofForGeometry(Geometry::Type GeomType) const;
496 
497  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
498  int Or) const;
499 
500  virtual const char * Name() const { return "QuadraticPos"; }
501 
502  virtual int GetContType() const { return CONTINUOUS; }
503 };
504 
505 /// Piecewise-(bi)cubic continuous finite elements.
507 {
508 private:
509  const PointFiniteElement PointFE;
510  const Cubic1DFiniteElement SegmentFE;
511  const Cubic2DFiniteElement TriangleFE;
512  const BiCubic2DFiniteElement QuadrilateralFE;
513  const Cubic3DFiniteElement TetrahedronFE;
514  const LagrangeHexFiniteElement ParallelepipedFE;
515  const H1_WedgeElement WedgeFE;
516 
517 public:
519  : ParallelepipedFE(3), WedgeFE(3, BasisType::ClosedUniform) { }
520 
521  virtual const FiniteElement *
523 
524  virtual int DofForGeometry(Geometry::Type GeomType) const;
525 
526  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
527  int Or) const;
528 
529  virtual const char * Name() const { return "Cubic"; }
530 
531  virtual int GetContType() const { return CONTINUOUS; }
532 };
533 
534 /// Crouzeix-Raviart nonconforming elements in 2D.
536 {
537 private:
538  const P0SegmentFiniteElement SegmentFE;
539  const CrouzeixRaviartFiniteElement TriangleFE;
540  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
541 public:
542  CrouzeixRaviartFECollection() : SegmentFE(1) { }
543 
544  virtual const FiniteElement *
546 
547  virtual int DofForGeometry(Geometry::Type GeomType) const;
548 
549  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
550  int Or) const;
551 
552  virtual const char * Name() const { return "CrouzeixRaviart"; }
553 
554  virtual int GetContType() const { return DISCONTINUOUS; }
555 };
556 
557 /// Piecewise-linear nonconforming finite elements in 3D.
559 {
560 private:
561  const P0TriangleFiniteElement TriangleFE;
562  const P1TetNonConfFiniteElement TetrahedronFE;
563  const P0QuadFiniteElement QuadrilateralFE;
564  const RotTriLinearHexFiniteElement ParallelepipedFE;
565 
566 public:
568 
569  virtual const FiniteElement *
571 
572  virtual int DofForGeometry(Geometry::Type GeomType) const;
573 
574  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
575  int Or) const;
576 
577  virtual const char * Name() const { return "LinearNonConf3D"; }
578 
579  virtual int GetContType() const { return DISCONTINUOUS; }
580 };
581 
582 
583 /** @brief First order Raviart-Thomas finite elements in 2D. This class is kept
584  only for backward compatibility, consider using RT_FECollection instead. */
586 {
587 private:
588  const P0SegmentFiniteElement SegmentFE; // normal component on edge
589  const RT0TriangleFiniteElement TriangleFE;
590  const RT0QuadFiniteElement QuadrilateralFE;
591 public:
592  RT0_2DFECollection() : SegmentFE(0) { }
593 
594  virtual const FiniteElement *
596 
597  virtual int DofForGeometry(Geometry::Type GeomType) const;
598 
599  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
600  int Or) const;
601 
602  virtual const char * Name() const { return "RT0_2D"; }
603 
604  virtual int GetContType() const { return NORMAL; }
605 };
606 
607 /** @brief Second order Raviart-Thomas finite elements in 2D. This class is kept
608  only for backward compatibility, consider using RT_FECollection instead. */
610 {
611 private:
612  const P1SegmentFiniteElement SegmentFE; // normal component on edge
613  const RT1TriangleFiniteElement TriangleFE;
614  const RT1QuadFiniteElement QuadrilateralFE;
615 public:
617 
618  virtual const FiniteElement *
620 
621  virtual int DofForGeometry(Geometry::Type GeomType) const;
622 
623  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
624  int Or) const;
625 
626  virtual const char * Name() const { return "RT1_2D"; }
627 
628  virtual int GetContType() const { return NORMAL; }
629 };
630 
631 /** @brief Third order Raviart-Thomas finite elements in 2D. This class is kept
632  only for backward compatibility, consider using RT_FECollection instead. */
634 {
635 private:
636  const P2SegmentFiniteElement SegmentFE; // normal component on edge
637  const RT2TriangleFiniteElement TriangleFE;
638  const RT2QuadFiniteElement QuadrilateralFE;
639 public:
641 
642  virtual const FiniteElement *
644 
645  virtual int DofForGeometry(Geometry::Type GeomType) const;
646 
647  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
648  int Or) const;
649 
650  virtual const char * Name() const { return "RT2_2D"; }
651 
652  virtual int GetContType() const { return NORMAL; }
653 };
654 
655 /** @brief Piecewise-constant discontinuous finite elements in 2D. This class is
656  kept only for backward compatibility, consider using L2_FECollection
657  instead. */
659 {
660 private:
661  const P0TriangleFiniteElement TriangleFE;
662  const P0QuadFiniteElement QuadrilateralFE;
663 public:
665 
666  virtual const FiniteElement *
668 
669  virtual int DofForGeometry(Geometry::Type GeomType) const;
670 
671  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
672  int Or) const;
673 
674  virtual const char * Name() const { return "Const2D"; }
675 
676  virtual int GetContType() const { return DISCONTINUOUS; }
677 };
678 
679 /** @brief Piecewise-linear discontinuous finite elements in 2D. This class is
680  kept only for backward compatibility, consider using L2_FECollection
681  instead. */
683 {
684 private:
685  const Linear2DFiniteElement TriangleFE;
686  const BiLinear2DFiniteElement QuadrilateralFE;
687 
688 public:
690 
691  virtual const FiniteElement *
693 
694  virtual int DofForGeometry(Geometry::Type GeomType) const;
695 
696  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
697  int Or) const;
698 
699  virtual const char * Name() const { return "LinearDiscont2D"; }
700 
701  virtual int GetContType() const { return DISCONTINUOUS; }
702 };
703 
704 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
706 {
707 private:
708  // const CrouzeixRaviartFiniteElement TriangleFE;
709  const GaussLinear2DFiniteElement TriangleFE;
710  const GaussBiLinear2DFiniteElement QuadrilateralFE;
711 
712 public:
714 
715  virtual const FiniteElement *
717 
718  virtual int DofForGeometry(Geometry::Type GeomType) const;
719 
720  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
721  int Or) const;
722 
723  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
724 
725  virtual int GetContType() const { return DISCONTINUOUS; }
726 };
727 
728 /// Linear (P1) finite elements on quadrilaterals.
730 {
731 private:
732  const P1OnQuadFiniteElement QuadrilateralFE;
733 public:
735  virtual const FiniteElement *
737  virtual int DofForGeometry(Geometry::Type GeomType) const;
738  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
739  int Or) const;
740  virtual const char * Name() const { return "P1OnQuad"; }
741  virtual int GetContType() const { return DISCONTINUOUS; }
742 };
743 
744 /** @brief Piecewise-quadratic discontinuous finite elements in 2D. This class
745  is kept only for backward compatibility, consider using L2_FECollection
746  instead. */
748 {
749 private:
750  const Quad2DFiniteElement TriangleFE;
751  const BiQuad2DFiniteElement QuadrilateralFE;
752 
753 public:
755 
756  virtual const FiniteElement *
758 
759  virtual int DofForGeometry(Geometry::Type GeomType) const;
760 
761  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
762  int Or) const;
763 
764  virtual const char * Name() const { return "QuadraticDiscont2D"; }
765  virtual int GetContType() const { return DISCONTINUOUS; }
766 };
767 
768 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
770 {
771 private:
772  const BiQuadPos2DFiniteElement QuadrilateralFE;
773 
774 public:
776  virtual const FiniteElement *
778  virtual int DofForGeometry(Geometry::Type GeomType) const;
779  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
780  int Or) const
781  { return NULL; }
782  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
783  virtual int GetContType() const { return DISCONTINUOUS; }
784 };
785 
786 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
788 {
789 private:
790  // const Quad2DFiniteElement TriangleFE;
791  const GaussQuad2DFiniteElement TriangleFE;
792  const GaussBiQuad2DFiniteElement QuadrilateralFE;
793 
794 public:
796 
797  virtual const FiniteElement *
799 
800  virtual int DofForGeometry(Geometry::Type GeomType) const;
801 
802  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
803  int Or) const;
804 
805  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
806  virtual int GetContType() const { return DISCONTINUOUS; }
807 };
808 
809 /** @brief Piecewise-cubic discontinuous finite elements in 2D. This class is
810  kept only for backward compatibility, consider using L2_FECollection
811  instead. */
813 {
814 private:
815  const Cubic2DFiniteElement TriangleFE;
816  const BiCubic2DFiniteElement QuadrilateralFE;
817 
818 public:
820 
821  virtual const FiniteElement *
823 
824  virtual int DofForGeometry(Geometry::Type GeomType) const;
825 
826  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
827  int Or) const;
828 
829  virtual const char * Name() const { return "CubicDiscont2D"; }
830  virtual int GetContType() const { return DISCONTINUOUS; }
831 };
832 
833 /** @brief Piecewise-constant discontinuous finite elements in 3D. This class is
834  kept only for backward compatibility, consider using L2_FECollection
835  instead. */
837 {
838 private:
839  const P0TetFiniteElement TetrahedronFE;
840  const P0HexFiniteElement ParallelepipedFE;
841  const L2_WedgeElement WedgeFE;
842 
843 public:
844  Const3DFECollection() : WedgeFE(0) { }
845 
846  virtual const FiniteElement *
848 
849  virtual int DofForGeometry(Geometry::Type GeomType) const;
850 
851  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
852  int Or) const;
853 
854  virtual const char * Name() const { return "Const3D"; }
855  virtual int GetContType() const { return DISCONTINUOUS; }
856 };
857 
858 /** @brief Piecewise-linear discontinuous finite elements in 3D. This class is
859  kept only for backward compatibility, consider using L2_FECollection
860  instead. */
862 {
863 private:
864  const Linear3DFiniteElement TetrahedronFE;
865  const TriLinear3DFiniteElement ParallelepipedFE;
866 
867 public:
869 
870  virtual const FiniteElement *
872 
873  virtual int DofForGeometry(Geometry::Type GeomType) const;
874 
875  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
876  int Or) const;
877 
878  virtual const char * Name() const { return "LinearDiscont3D"; }
879  virtual int GetContType() const { return DISCONTINUOUS; }
880 };
881 
882 /** @brief Piecewise-quadratic discontinuous finite elements in 3D. This class
883  is kept only for backward compatibility, consider using L2_FECollection
884  instead. */
886 {
887 private:
888  const Quadratic3DFiniteElement TetrahedronFE;
889  const LagrangeHexFiniteElement ParallelepipedFE;
890 
891 public:
892  QuadraticDiscont3DFECollection () : ParallelepipedFE(2) { }
893 
894  virtual const FiniteElement *
896 
897  virtual int DofForGeometry(Geometry::Type GeomType) const;
898 
899  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
900  int Or) const;
901 
902  virtual const char * Name() const { return "QuadraticDiscont3D"; }
903  virtual int GetContType() const { return DISCONTINUOUS; }
904 };
905 
906 /// Finite element collection on a macro-element.
908 {
909 private:
910  const PointFiniteElement PointFE;
911  const RefinedLinear1DFiniteElement SegmentFE;
912  const RefinedLinear2DFiniteElement TriangleFE;
913  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
914  const RefinedLinear3DFiniteElement TetrahedronFE;
915  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
916 
917 public:
919 
920  virtual const FiniteElement *
922 
923  virtual int DofForGeometry(Geometry::Type GeomType) const;
924 
925  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
926  int Or) const;
927 
928  virtual const char * Name() const { return "RefinedLinear"; }
929  virtual int GetContType() const { return CONTINUOUS; }
930 };
931 
932 /** @brief Lowest order Nedelec finite elements in 3D. This class is kept only
933  for backward compatibility, consider using the new ND_FECollection
934  instead. */
936 {
937 private:
938  const Nedelec1HexFiniteElement HexahedronFE;
939  const Nedelec1TetFiniteElement TetrahedronFE;
940 
941 public:
943 
944  virtual const FiniteElement *
946 
947  virtual int DofForGeometry(Geometry::Type GeomType) const;
948 
949  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
950  int Or) const;
951 
952  virtual const char * Name() const { return "ND1_3D"; }
953  virtual int GetContType() const { return TANGENTIAL; }
954 };
955 
956 /** @brief First order Raviart-Thomas finite elements in 3D. This class is kept
957  only for backward compatibility, consider using RT_FECollection instead. */
959 {
960 private:
961  const P0TriangleFiniteElement TriangleFE;
962  const P0QuadFiniteElement QuadrilateralFE;
963  const RT0HexFiniteElement HexahedronFE;
964  const RT0TetFiniteElement TetrahedronFE;
965 public:
967 
968  virtual const FiniteElement *
970 
971  virtual int DofForGeometry(Geometry::Type GeomType) const;
972 
973  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
974  int Or) const;
975 
976  virtual const char * Name() const { return "RT0_3D"; }
977  virtual int GetContType() const { return NORMAL; }
978 };
979 
980 /** @brief Second order Raviart-Thomas finite elements in 3D. This class is kept
981  only for backward compatibility, consider using RT_FECollection instead. */
983 {
984 private:
985  const Linear2DFiniteElement TriangleFE;
986  const BiLinear2DFiniteElement QuadrilateralFE;
987  const RT1HexFiniteElement HexahedronFE;
988 public:
990 
991  virtual const FiniteElement *
993 
994  virtual int DofForGeometry(Geometry::Type GeomType) const;
995 
996  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
997  int Or) const;
998 
999  virtual const char * Name() const { return "RT1_3D"; }
1000  virtual int GetContType() const { return NORMAL; }
1001 };
1002 
1003 /// Discontinuous collection defined locally by a given finite element.
1005 {
1006 private:
1007  char d_name[32];
1008  Geometry::Type GeomType;
1009  FiniteElement *Local_Element;
1010 
1011 public:
1012  Local_FECollection(const char *fe_name);
1013 
1015  Geometry::Type _GeomType) const
1016  { return (GeomType == _GeomType) ? Local_Element : NULL; }
1017  virtual int DofForGeometry(Geometry::Type _GeomType) const
1018  { return (GeomType == _GeomType) ? Local_Element->GetDof() : 0; }
1019  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1020  int Or) const
1021  { return NULL; }
1022  virtual const char *Name() const { return d_name; }
1023 
1024  virtual ~Local_FECollection() { delete Local_Element; }
1025  virtual int GetContType() const { return DISCONTINUOUS; }
1026 };
1027 
1028 }
1029 
1030 #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 all finite elements.
Definition: fe.hpp:235
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:403
Normal component of vector field.
Definition: fe_coll.hpp:47
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
A 2D 2nd order Raviart-Thomas vector element on a triangle.
Definition: fe.hpp:1396
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:372
Field is discontinuous across element interfaces.
Definition: fe_coll.hpp:48
virtual const char * Name() const
Definition: fe_coll.hpp:529
virtual const char * Name() const
Definition: fe_coll.hpp:902
virtual const char * Name() const
Definition: fe_coll.hpp:805
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1887
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1169
Arbitrary order L2 elements in 3D on a wedge.
Definition: fe.hpp:2619
Piecewise-cubic discontinuous finite elements in 2D. This class is kept only for backward compatibili...
Definition: fe_coll.hpp:812
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
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:787
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:278
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:337
A 2D 1st order Raviart-Thomas vector element on a square.
Definition: fe.hpp:1367
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:779
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1210
A 1D constant element on a segment.
Definition: fe.hpp:1327
A 1D linear element with nodes at 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:1505
A 1D quadratic element with nodes at the Gaussian points (trace of RT2)
Definition: fe.hpp:1516
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:928
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:854
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:173
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:2777
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:2167
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:2262
A 3D 1st order Nedelec element on a tetrahedron.
Definition: fe.hpp:1721
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:699
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:2393
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2811
static const int NumGeom
Definition: geom.hpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:477
virtual int GetContType() const =0
Tangential components of vector field.
Definition: fe_coll.hpp:46
const Geometry::Type geom
Definition: ex1.cpp:40
A 2D 3rd order Raviart-Thomas vector element on a square.
Definition: fe.hpp:1474
virtual ~L2_FECollection()
Definition: fe_coll.cpp:2186
virtual const char * Name() const
Definition: fe_coll.hpp:740
virtual int GetContType() const
Definition: fe_coll.hpp:479
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2785
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/tri)linear continuous finite elements.
Definition: fe_coll.hpp:428
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:1812
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:928
void Reset() const
Definition: fe.hpp:3215
A 2D quadratic element on triangle with nodes at the vertices and midpoints of the triangle...
Definition: fe.hpp:1048
A 0D point finite element.
Definition: fe.hpp:896
A 3D 0th order Raviert-Thomas element on a cube.
Definition: fe.hpp:1745
A 1D quadractic finite element with uniformly spaced nodes.
Definition: fe.hpp:1016
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1377
virtual const char * Name() const
Definition: fe_coll.hpp:303
virtual int GetContType() const
Definition: fe_coll.hpp:929
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2661
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
Second order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:609
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:176
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:2275
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2411
virtual const char * Name() const
Definition: fe_coll.hpp:417
A 3D constant element on a tetrahedron.
Definition: fe.hpp:1554
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:2459
int HasFaceDofs(Geometry::Type GeomType) const
Definition: fe_coll.cpp:25
virtual const char * Name() const
Definition: fe_coll.hpp:602
A 3D Crouzeix-Raviart element on the tetrahedron.
Definition: fe.hpp:1543
A 2D bi-cubic element on a square with uniformly spaces nodes.
Definition: fe.hpp:1145
virtual const char * Name() const
Definition: fe_coll.hpp:829
Finite element collection on a macro-element.
Definition: fe_coll.hpp:907
virtual int GetContType() const
Definition: fe_coll.hpp:903
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
A 2D refined bi-linear FE on a square.
Definition: fe.hpp:1656
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:2754
virtual int GetContType() const
Definition: fe_coll.hpp:181
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1184
virtual int GetContType() const
Definition: fe_coll.hpp:855
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:353
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:2427
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1907
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:238
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:506
Class for finite elements with basis functions that return scalar values.
Definition: fe.hpp:615
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:824
virtual const char * Name() const
Definition: fe_coll.hpp:626
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1853
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
virtual int GetContType() const
Definition: fe_coll.hpp:879
virtual int GetContType() const
Definition: fe_coll.hpp:304
virtual int GetContType() const
Definition: fe_coll.hpp:652
A 2D bi-linear element on a square with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:987
A 1D refined linear element.
Definition: fe.hpp:1603
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:166
A 2D refined linear element on a triangle.
Definition: fe.hpp:1623
A 2D refined linear element on a tetrahedron.
Definition: fe.hpp:1643
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:277
virtual int GetContType() const
Definition: fe_coll.hpp:701
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:906
A 2D constant element on a triangle.
Definition: fe.hpp:1203
virtual int DofForGeometry(Geometry::Type _GeomType) const
Definition: fe_coll.hpp:1017
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:293
int GetBasisType() const
Definition: fe_coll.hpp:263
virtual int GetContType() const
Definition: fe_coll.hpp:953
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1365
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
Piecewise-linear discontinuous finite elements in 2D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:682
int GetBasisType() const
Definition: fe_coll.hpp:184
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:2693
Arbitrary order &quot;H^{1/2}-conforming&quot; trace finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:213
virtual int GetContType() const
Definition: fe_coll.hpp:354
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2479
A 3D refined tri-linear element on a cube.
Definition: fe.hpp:1676
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:1004
Arbitrary order H(curl)-trace finite elements defined on the interface between mesh elements (faces...
Definition: fe_coll.hpp:363
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:769
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2681
A 3D quadratic element on a tetrahedron with uniformly spaced nodes.
Definition: fe.hpp:1265
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:811
virtual int GetContType() const
Definition: fe_coll.hpp:628
Third order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:633
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:963
First order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:958
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:483
virtual ~Local_FECollection()
Definition: fe_coll.hpp:1024
virtual const char * Name() const
Definition: fe_coll.hpp:723
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:67
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:582
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1874
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
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:836
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:1019
A 3D linear element on a tetrahedron with nodes at the vertices of the tetrahedron.
Definition: fe.hpp:1236
A 1D linear element with nodes on the endpoints.
Definition: fe.hpp:909
virtual int GetContType() const
Definition: fe_coll.hpp:255
A 2D Crouzeix-Raviart finite element on square.
Definition: fe.hpp:1315
A 2D linear element on triangle with nodes at the vertices of the triangle.
Definition: fe.hpp:929
virtual int GetContType() const
Definition: fe_coll.hpp:765
A quadratic element on triangle with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:1072
virtual const char * Name() const
Definition: fe_coll.hpp:180
virtual const char * Name() const
Definition: fe_coll.hpp:878
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1834
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:299
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:658
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:165
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:558
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:535
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:272
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:269
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:196
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:243
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:257
H1Ser_FECollection(const int p, const int dim=2)
Definition: fe_coll.hpp:206
A 2D 2nd order Raviart-Thomas vector element on a square.
Definition: fe.hpp:1425
virtual const char * Name() const
Definition: fe_coll.hpp:552
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:449
virtual const char * Name() const
Definition: fe_coll.hpp:650
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:2764
virtual const char * Name() const
Definition: fe_coll.hpp:577
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:865
A 3D 1st order Raviert-Thomas element on a cube.
Definition: fe.hpp:1775
virtual const char * Name() const
Definition: fe_coll.hpp:952
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
A 2D bi-quadratic element on a square with uniformly spaced nodes.
Definition: fe.hpp:1089
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:2638
A 3D tri-linear element on a cube with nodes at the vertices of the cube.
Definition: fe.hpp:1278
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type _GeomType) const
Definition: fe_coll.hpp:1014
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:745
Piecewise-quadratic discontinuous finite elements in 3D. This class is kept only for backward compati...
Definition: fe_coll.hpp:885
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:599
Field is continuous across element interfaces.
Definition: fe_coll.hpp:45
virtual int GetContType() const
Definition: fe_coll.hpp:977
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
A 3D constant element on a cube.
Definition: fe.hpp:1567
virtual const char * Name() const
Definition: fe_coll.hpp:764
A 2D bi-linear element on a square with nodes at the vertices of the square.
Definition: fe.hpp:951
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:705
virtual const char * Name() const
Definition: fe_coll.hpp:253
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:47
A linear element on a triangle with nodes at the 3 &quot;Gaussian&quot; points.
Definition: fe.hpp:975
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1322
virtual int GetContType() const
Definition: fe_coll.hpp:554
virtual int GetContType() const
Definition: fe_coll.hpp:531
Lowest order Nedelec finite elements in 3D. This class is kept only for backward compatibility, consider using the new ND_FECollection instead.
Definition: fe_coll.hpp:935
Arbitrary order H1-conforming (continuous) finite elements with positive basis functions.
Definition: fe_coll.hpp:193
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:61
virtual int GetContType() const
Definition: fe_coll.hpp:676
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1469
virtual int GetContType() const
Definition: fe_coll.hpp:502
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:349
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
A 2D bi-quadratic element on a square with nodes at the 9 &quot;Gaussian&quot; points.
Definition: fe.hpp:1132
virtual const char * Name() const
Definition: fe_coll.hpp:674
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:2804
A 1D quadratic positive element utilizing the 2nd order Bernstein basis.
Definition: fe.hpp:1036
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:639
int dim
Definition: ex24.cpp:53
Piecewise-linear discontinuous finite elements in 3D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:861
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:296
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:976
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:679
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1422
virtual int GetContType() const
Definition: fe_coll.hpp:419
virtual int GetContType() const
Definition: fe_coll.hpp:725
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:2439
virtual int GetContType() const
Definition: fe_coll.hpp:830
Tensor products of 1D Lagrange1DFiniteElement (only degree 2 is functional)
Definition: fe.hpp:1581
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:557
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:455
First order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:585
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
A 2D constant element on a square.
Definition: fe.hpp:1221
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:333
void Reset() const
Definition: fe_coll.hpp:394
A 2D cubic element on a triangle with uniformly spaced nodes.
Definition: fe.hpp:1173
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1044
A 1D cubic element with uniformly spaced nodes.
Definition: fe.hpp:1160
virtual const char * Name() const
Definition: fe_coll.hpp:1022
Arbitrary order &quot;H^{-1/2}-conforming&quot; face finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:313
virtual int GetContType() const
Definition: fe_coll.hpp:1025
A 2D 1st order Raviart-Thomas vector element on a triangle.
Definition: fe.hpp:1338
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
Second order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:982
virtual const char * Name() const
Definition: fe_coll.hpp:999
A 2D linear element on a square with 3 nodes at the vertices of the lower left triangle.
Definition: fe.hpp:1003
virtual const char * Name() const
Definition: fe_coll.hpp:782
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1095
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2798
Piecewise-quadratic discontinuous finite elements in 2D. This class is kept only for backward compati...
Definition: fe_coll.hpp:747
Arbitrary order H1 elements in 3D on a wedge.
Definition: fe.hpp:2356
A 3D 0th order Raviert-Thomas element on a tetrahedron.
Definition: fe.hpp:1805
virtual int GetContType() const
Definition: fe_coll.hpp:579
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
A 2D 3rd order Raviart-Thomas vector element on a triangle.
Definition: fe.hpp:1454
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:2712
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1147
virtual const char * Name() const
Definition: fe_coll.hpp:500
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:729
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
An arbitrary order 3D NURBS element on a cube.
Definition: fe.hpp:3281
virtual int GetContType() const
Definition: fe_coll.hpp:741
An arbitrary order 1D NURBS element on a segment.
Definition: fe.hpp:3229
virtual int GetContType() const
Definition: fe_coll.hpp:1000
virtual ~FiniteElementCollection()
Definition: fe_coll.hpp:75
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1134
virtual int GetContType() const
Definition: fe_coll.hpp:604
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:338
An arbitrary order 2D NURBS element on a square.
Definition: fe.hpp:3249
A 3D 1st order Nedelec element on a cube.
Definition: fe.hpp:1697
virtual int GetContType() const
Definition: fe_coll.hpp:451
A 2D Crouzeix-Raviart element on triangle.
Definition: fe.hpp:1302
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:221
double f(const Vector &p)
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:346