MFEM  v4.6.0
Finite element discretization library
fe_coll.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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  /** @brief Returns the first non-NULL FiniteElement for the given dimension
55 
56  @note Repeatedly calls FiniteElementForGeometry in the order defined in
57  the Geometry::Type enumeration.
58  */
59  virtual const FiniteElement *
60  FiniteElementForDim(int dim) const;
61 
62  virtual int DofForGeometry(Geometry::Type GeomType) const = 0;
63 
64  /** @brief Returns a DoF transformation object compatible with this basis
65  and geometry type.
66  */
69  { return NULL; }
70 
71  /** @brief Returns an array, say p, that maps a local permuted index i to a
72  local base index: base_i = p[i].
73 
74  @note Only provides information about interior dofs. See
75  FiniteElementCollection::SubDofOrder if interior \a and boundary dof
76  order is needed. */
77  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
78  int Or) const = 0;
79 
80  virtual const char * Name() const { return "Undefined"; }
81 
82  virtual int GetContType() const = 0;
83 
84  /** @note The following methods provide the same information as the
85  corresponding methods of the FiniteElement base class.
86  @{
87  */
88  virtual int GetRangeType(int dim) const;
89  virtual int GetDerivRangeType(int dim) const;
90  virtual int GetMapType(int dim) const;
91  virtual int GetDerivType(int dim) const;
92  virtual int GetDerivMapType(int dim) const;
93  /** @} */
94 
95  int HasFaceDofs(Geometry::Type geom, int p) const;
96 
98  Geometry::Type GeomType) const
99  {
100  return FiniteElementForGeometry(GeomType);
101  }
102 
104 
105  virtual ~FiniteElementCollection();
106 
107  /** @brief Factory method: return a newly allocated FiniteElementCollection
108  according to the given name. */
109  /**
110  | FEC Name | Space | Order | BasisType | FiniteElement::MapT | Notes |
111  | :------: | :---: | :---: | :-------: | :-----: | :---: |
112  | H1_[DIM]_[ORDER] | H1 | * | 1 | VALUE | H1 nodal elements |
113  | H1@[BTYPE]_[DIM]_[ORDER] | H1 | * | * | VALUE | H1 nodal elements |
114  | H1Pos_[DIM]_[ORDER] | H1 | * | 1 | VALUE | H1 nodal elements |
115  | 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) |
116  | 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) |
117  | 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) |
118  | ND_[DIM]_[ORDER] | H(curl) | * | 1 / 0 | H_CURL | Nedelec vector elements |
119  | ND@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(curl) | * | * / * | H_CURL | Nedelec vector elements |
120  | 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) |
121  | 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) |
122  | RT_[DIM]_[ORDER] | H(div) | * | 1 / 0 | H_DIV | Raviart-Thomas vector elements |
123  | RT@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(div) | * | * / * | H_DIV | Raviart-Thomas vector elements |
124  | 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) |
125  | 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) |
126  | 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) |
127  | 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) |
128  | L2_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinuous L2 elements |
129  | L2_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinuous L2 elements |
130  | L2Int_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinuous L2 elements |
131  | L2Int_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinuous L2 elements |
132  | DG_Iface_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
133  | DG_Iface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
134  | DG_IntIface_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
135  | DG_IntIface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
136  | NURBS[ORDER] | - | * | - | VALUE | Non-Uniform Rational B-Splines (NURBS) elements |
137  | LinearNonConf3D | - | 1 | 1 | VALUE | Piecewise-linear nonconforming finite elements in 3D |
138  | CrouzeixRaviart | - | - | - | - | Crouzeix-Raviart nonconforming elements in 2D |
139  | Local_[FENAME] | - | - | - | - | Special collection that builds a local version out of the FENAME collection |
140  |-|-|-|-|-|-|
141  | Linear | H1 | 1 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
142  | Quadratic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
143  | QuadraticPos | H1 | 2 | 2 | VALUE | Left in for backward compatibility, consider using H1_ |
144  | Cubic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
145  | Const2D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
146  | Const3D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
147  | LinearDiscont2D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
148  | GaussLinearDiscont2D | L2 | 1 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
149  | P1OnQuad | H1 | 1 | 1 | VALUE | Linear P1 element with 3 nodes on a square |
150  | QuadraticDiscont2D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
151  | QuadraticPosDiscont2D | L2 | 2 | 2 | VALUE | Left in for backward compatibility, consider using L2_ |
152  | GaussQuadraticDiscont2D | L2 | 2 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
153  | CubicDiscont2D | L2 | 3 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
154  | LinearDiscont3D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
155  | QuadraticDiscont3D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
156  | ND1_3D | H(Curl) | 1 | 1 / 0 | H_CURL | Left in for backward compatibility, consider using ND_ |
157  | RT0_2D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
158  | RT1_2D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
159  | RT2_2D | H(Div) | 3 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
160  | RT0_3D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
161  | RT1_3D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
162 
163  | Tag | Description |
164  | :------: | :--------: |
165  | [DIM] | Dimension of the elements (1D, 2D, 3D) |
166  | [ORDER] | Approximation order of the elements (P0, P1, P2, ...) |
167  | [BTYPE] | BasisType of the element (0-GaussLegendre, 1 - GaussLobatto, 2-Bernstein, 3-OpenUniform, 4-CloseUniform, 5-OpenHalfUniform) |
168  | [OBTYPE] | Open BasisType of the element for elements which have both types |
169  | [CBTYPE] | Closed BasisType of the element for elements which have both types |
170 
171  [FENAME] Is a special case for the Local FEC which generates a local version of a given
172  FEC. It is selected from one of (BiCubic2DFiniteElement, Quad_Q3, Nedelec1HexFiniteElement,
173  Hex_ND1, H1_[DIM]_[ORDER],H1Pos_[DIM]_[ORDER], L2_[DIM]_[ORDER] )
174  */
175  static FiniteElementCollection *New(const char *name);
176 
177  /** @brief Get the local dofs for a given sub-manifold.
178 
179  Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex, 1D
180  - edge, 2D - face) including those on its boundary. The local index of the
181  sub-manifold (inside Geom) and its orientation are given by the parameter
182  Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed that 0 <=
183  SDim <= Dim(Geom). */
184  void SubDofOrder(Geometry::Type Geom, int SDim, int Info,
185  Array<int> &dofs) const;
186 
187  /// Variable order version of FiniteElementForGeometry().
188  /** The order parameter @a p represents the order of the highest-dimensional
189  FiniteElement%s the fixed-order collection we want to query. In general,
190  this order is different from the order of the returned FiniteElement. */
191  const FiniteElement *GetFE(Geometry::Type geom, int p) const
192  {
193  if (p == base_p) { return FiniteElementForGeometry(geom); }
194  if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
195  return var_orders[p]->FiniteElementForGeometry(geom);
196  }
197 
198  /// Variable order version of DofForGeometry().
199  /** The order parameter @a p represents the order of the highest-dimensional
200  FiniteElement%s the fixed-order collection we want to query. In general,
201  this order is different from the order of the element corresponding to
202  @a geom in that fixed-order collection. */
203  int GetNumDof(Geometry::Type geom, int p) const
204  {
205  if (p == base_p) { return DofForGeometry(geom); }
206  if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
207  return var_orders[p]->DofForGeometry(geom);
208  }
209 
210  /// Variable order version of DofOrderForOrientation().
211  /** The order parameter @a p represents the order of the highest-dimensional
212  FiniteElement%s the fixed-order collection we want to query. In general,
213  this order is different from the order of the element corresponding to
214  @a geom in that fixed-order collection. */
215  const int *GetDofOrdering(Geometry::Type geom, int p, int ori) const
216  {
217  if (p == base_p) { return DofOrderForOrientation(geom, ori); }
218  if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
219  return var_orders[p]->DofOrderForOrientation(geom, ori);
220  }
221 
222  /** @brief Return the order (polynomial degree) of the FE collection,
223  corresponding to the order/degree returned by FiniteElement::GetOrder()
224  of the highest-dimensional FiniteElement%s defined by the collection. */
225  int GetOrder() const { return base_p; }
226 
227  /// Instantiate a new collection of the same type with a different order.
228  /** Generally, the order parameter @a p is NOT the same as the parameter @a p
229  used by some of the constructors of derived classes. Instead, this @a p
230  represents the order of the new FE collection as it will be returned by
231  its GetOrder() method. */
232  virtual FiniteElementCollection *Clone(int p) const;
233 
234 protected:
235  const int base_p; ///< Order as returned by GetOrder().
236 
239 
240  void InitVarOrder(int p) const;
241 
243 
244  /// How to treat errors in FiniteElementForGeometry() calls.
246  {
247  RETURN_NULL, ///< Return NULL on errors
248  RAISE_MFEM_ERROR /**< Raise an MFEM error (default in base class).
249  Sub-classes can ignore this and return NULL. */
250  };
251 
252  /// How to treat errors in FiniteElementForGeometry() calls.
253  /** The typical error in derived classes is that no FiniteElement is defined
254  for the given Geometry, or the input is not a valid Geometry. */
256 };
257 
258 /// Arbitrary order H1-conforming (continuous) finite elements.
260 {
261 
262 protected:
263  int dim, b_type;
264  char h1_name[32];
267  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8], *TetDofOrd[24];
268 
269 public:
270  explicit H1_FECollection(const int p, const int dim = 3,
271  const int btype = BasisType::GaussLobatto);
272 
274  Geometry::Type GeomType) const;
275  virtual int DofForGeometry(Geometry::Type GeomType) const
276  { return H1_dof[GeomType]; }
277  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
278  int Or) const;
279 
280  virtual const char *Name() const { return h1_name; }
281  virtual int GetContType() const { return CONTINUOUS; }
282  int GetBasisType() const { return b_type; }
283 
285 
286  /// Get the Cartesian to local H1 dof map
287  const int *GetDofMap(Geometry::Type GeomType) const;
288  /// Variable order version of GetDofMap
289  const int *GetDofMap(Geometry::Type GeomType, int p) const;
290 
292  { return new H1_FECollection(p, dim, b_type); }
293 
294  virtual ~H1_FECollection();
295 };
296 
297 /** @brief Arbitrary order H1-conforming (continuous) finite elements with
298  positive basis functions. */
300 {
301 public:
302  explicit H1Pos_FECollection(const int p, const int dim = 3)
303  : H1_FECollection(p, dim, BasisType::Positive) { }
304 };
305 
306 
307 /** Arbitrary order H1-conforming (continuous) serendipity finite elements;
308  Current implementation works in 2D only; 3D version is in development. */
310 {
311 public:
312  explicit H1Ser_FECollection(const int p, const int dim = 2)
313  : H1_FECollection(p, dim, BasisType::Serendipity) { };
314 };
315 
316 /** @brief Arbitrary order "H^{1/2}-conforming" trace finite elements defined on
317  the interface between mesh elements (faces,edges,vertices); these are the
318  trace FEs of the H1-conforming FEs. */
320 {
321 public:
322  H1_Trace_FECollection(const int p, const int dim,
323  const int btype = BasisType::GaussLobatto);
324 };
325 
326 /// Arbitrary order "L2-conforming" discontinuous finite elements.
328 {
329 private:
330  int dim;
331  int b_type; // BasisType
332  int m_type; // map type
333  char d_name[32];
336  int *SegDofOrd[2]; // for rotating segment dofs in 1D
337  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
338  int *TetDofOrd[24]; // for rotating tetrahedron dofs in 3D
339  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
340 
341 public:
342  L2_FECollection(const int p, const int dim,
343  const int btype = BasisType::GaussLegendre,
344  const int map_type = FiniteElement::VALUE);
345 
347  Geometry::Type GeomType) const;
348  virtual int DofForGeometry(Geometry::Type GeomType) const
349  {
350  if (L2_Elements[GeomType])
351  {
352  return L2_Elements[GeomType]->GetDof();
353  }
354  return 0;
355  }
356  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
357  int Or) const;
358  virtual const char *Name() const { return d_name; }
359 
360  virtual int GetContType() const { return DISCONTINUOUS; }
361 
363  Geometry::Type GeomType) const
364  {
365  return Tr_Elements[GeomType];
366  }
367 
368  int GetBasisType() const { return b_type; }
369 
371  { return new L2_FECollection(p, dim, b_type, m_type); }
372 
373  virtual ~L2_FECollection();
374 };
375 
376 /// Declare an alternative name for L2_FECollection = DG_FECollection
378 
379 /// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
381 {
382 protected:
383  int dim;
384  int cb_type; // closed BasisType
385  int ob_type; // open BasisType
386  char rt_name[32];
389  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
390 
391  // Initialize only the face elements
392  void InitFaces(const int p, const int dim, const int map_type,
393  const bool signs);
394 
395  // Constructor used by the constructor of the RT_Trace_FECollection and
396  // DG_Interface_FECollection classes
397  RT_FECollection(const int p, const int dim, const int map_type,
398  const bool signs,
399  const int ob_type = BasisType::GaussLegendre);
400 
401 public:
402  /// Construct an H(div)-conforming Raviart-Thomas FE collection, RT_p.
403  /** The index @a p corresponds to the space RT_p, as typically denoted in the
404  literature, which contains vector polynomials of degree up to (p+1).
405  For example, the RT_0 collection contains vector-valued linear functions
406  and, in particular, FiniteElementCollection::GetOrder() will,
407  correspondingly, return order 1. */
408  RT_FECollection(const int p, const int dim,
409  const int cb_type = BasisType::GaussLobatto,
410  const int ob_type = BasisType::GaussLegendre);
411 
413  Geometry::Type GeomType) const;
414  virtual int DofForGeometry(Geometry::Type GeomType) const
415  { return RT_dof[GeomType]; }
416  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
417  int Or) const;
418  virtual const char *Name() const { return rt_name; }
419  virtual int GetContType() const { return NORMAL; }
421 
422  int GetClosedBasisType() const { return cb_type; }
423  int GetOpenBasisType() const { return ob_type; }
424 
426  { return new RT_FECollection(p, dim, cb_type, ob_type); }
427 
428  virtual ~RT_FECollection();
429 };
430 
431 /** @brief Arbitrary order "H^{-1/2}-conforming" face finite elements defined on
432  the interface between mesh elements (faces); these are the normal trace FEs
433  of the H(div)-conforming FEs. */
435 {
436 public:
437  RT_Trace_FECollection(const int p, const int dim,
438  const int map_type = FiniteElement::INTEGRAL,
439  const int ob_type = BasisType::GaussLegendre);
440 };
441 
442 /** Arbitrary order discontinuous finite elements defined on the interface
443  between mesh elements (faces). The functions in this space are single-valued
444  on each face and are discontinuous across its boundary. */
446 {
447 public:
448  DG_Interface_FECollection(const int p, const int dim,
449  const int map_type = FiniteElement::VALUE,
450  const int ob_type = BasisType::GaussLegendre);
451 };
452 
453 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
455 {
456 protected:
457  int dim;
458  int cb_type; // closed BasisType
459  int ob_type; // open BasisType
460  char nd_name[32];
463  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
464 
465 public:
466  ND_FECollection(const int p, const int dim,
467  const int cb_type = BasisType::GaussLobatto,
468  const int ob_type = BasisType::GaussLegendre);
469 
470  virtual const FiniteElement *
472 
473  virtual int DofForGeometry(Geometry::Type GeomType) const
474  { return ND_dof[GeomType]; }
475 
478 
479  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
480  int Or) const;
481 
482  virtual const char *Name() const { return nd_name; }
483  virtual int GetContType() const { return TANGENTIAL; }
485 
486  int GetClosedBasisType() const { return cb_type; }
487  int GetOpenBasisType() const { return ob_type; }
488 
490  { return new ND_FECollection(p, dim, cb_type, ob_type); }
491 
492  virtual ~ND_FECollection();
493 };
494 
495 /** @brief Arbitrary order H(curl)-trace finite elements defined on the
496  interface between mesh elements (faces,edges); these are the tangential
497  trace FEs of the H(curl)-conforming FEs. */
499 {
500 public:
501  ND_Trace_FECollection(const int p, const int dim,
502  const int cb_type = BasisType::GaussLobatto,
503  const int ob_type = BasisType::GaussLegendre);
504 };
505 
506 /// Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
508 {
509 protected:
510  char nd_name[32];
513 
514 public:
515  ND_R1D_FECollection(const int p, const int dim,
516  const int cb_type = BasisType::GaussLobatto,
517  const int ob_type = BasisType::GaussLegendre);
518 
520  const
521  { return ND_Elements[GeomType]; }
522  virtual int DofForGeometry(Geometry::Type GeomType) const
523  { return ND_dof[GeomType]; }
524  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
525  int Or) const;
526  virtual const char *Name() const { return nd_name; }
527  virtual int GetContType() const { return TANGENTIAL; }
529 
530  virtual ~ND_R1D_FECollection();
531 };
532 
533 /// Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 1D.
535 {
536 protected:
537  char rt_name[32];
540 
541 public:
542  RT_R1D_FECollection(const int p, const int dim,
543  const int cb_type = BasisType::GaussLobatto,
544  const int ob_type = BasisType::GaussLegendre);
545 
547  const
548  { return RT_Elements[GeomType]; }
549  virtual int DofForGeometry(Geometry::Type GeomType) const
550  { return RT_dof[GeomType]; }
551  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
552  int Or) const;
553  virtual const char *Name() const { return rt_name; }
554  virtual int GetContType() const { return NORMAL; }
556 
557  virtual ~RT_R1D_FECollection();
558 };
559 
560 /// Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
562 {
563 protected:
564  char nd_name[32];
567  int *SegDofOrd[2];
568 
569 public:
570  ND_R2D_FECollection(const int p, const int dim,
571  const int cb_type = BasisType::GaussLobatto,
572  const int ob_type = BasisType::GaussLegendre);
573 
575  const
576  { return ND_Elements[GeomType]; }
577  virtual int DofForGeometry(Geometry::Type GeomType) const
578  { return ND_dof[GeomType]; }
579  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
580  int Or) const;
581  virtual const char *Name() const { return nd_name; }
582  virtual int GetContType() const { return TANGENTIAL; }
584 
585  virtual ~ND_R2D_FECollection();
586 };
587 
588 /** @brief Arbitrary order 3D H(curl)-trace finite elements in 2D defined on the
589  interface between mesh elements (edges); these are the tangential
590  trace FEs of the H(curl)-conforming FEs. */
592 {
593 public:
594  ND_R2D_Trace_FECollection(const int p, const int dim,
595  const int cb_type = BasisType::GaussLobatto,
596  const int ob_type = BasisType::GaussLegendre);
597 };
598 
599 /// Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 2D.
601 {
602 protected:
603  int ob_type; // open BasisType
604  char rt_name[32];
607  int *SegDofOrd[2];
608 
609  // Initialize only the face elements
610  void InitFaces(const int p, const int dim, const int map_type,
611  const bool signs);
612 
613  // Constructor used by the constructor of the RT_R2D_Trace_FECollection
614  RT_R2D_FECollection(const int p, const int dim, const int map_type,
615  const bool signs,
616  const int ob_type = BasisType::GaussLegendre);
617 
618 public:
619  RT_R2D_FECollection(const int p, const int dim,
620  const int cb_type = BasisType::GaussLobatto,
621  const int ob_type = BasisType::GaussLegendre);
622 
624  const
625  { return RT_Elements[GeomType]; }
626  virtual int DofForGeometry(Geometry::Type GeomType) const
627  { return RT_dof[GeomType]; }
628  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
629  int Or) const;
630  virtual const char *Name() const { return rt_name; }
631  virtual int GetContType() const { return NORMAL; }
633 
634  virtual ~RT_R2D_FECollection();
635 };
636 
637 /** @brief Arbitrary order 3D "H^{-1/2}-conforming" face finite elements defined on
638  the interface between mesh elements (faces); these are the normal trace FEs
639  of the H(div)-conforming FEs. */
641 {
642 public:
643  RT_R2D_Trace_FECollection(const int p, const int dim,
644  const int map_type = FiniteElement::INTEGRAL,
645  const int ob_type = BasisType::GaussLegendre);
646 };
647 
648 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
650 {
651 private:
652  PointFiniteElement *PointFE;
653  NURBS1DFiniteElement *SegmentFE;
654  NURBS2DFiniteElement *QuadrilateralFE;
655  NURBS3DFiniteElement *ParallelepipedFE;
656 
657  mutable int mOrder; // >= 1 or VariableOrder
658  // The 'name' can be:
659  // 1) name = "NURBS" + "number", for fixed order, or
660  // 2) name = "NURBS", for VariableOrder.
661  // The name is updated before writing it to a stream, for example, see
662  // FiniteElementSpace::Save().
663  mutable char name[16];
664 
665 public:
666  enum { VariableOrder = -1 };
667 
668  /** @brief The parameter @a Order must be either a positive number, for fixed
669  order, or VariableOrder (default). */
670  explicit NURBSFECollection(int Order = VariableOrder);
671 
672  void Reset() const
673  {
674  SegmentFE->Reset();
675  QuadrilateralFE->Reset();
676  ParallelepipedFE->Reset();
677  }
678 
679  /** @brief Get the order of the NURBS collection: either a positive number,
680  when using fixed order, or VariableOrder. */
681  /** @note Not to be confused with FiniteElementCollection::GetOrder(). */
682  int GetOrder() const { return mOrder; }
683 
684  /** @brief Set the order and the name, based on the given @a Order: either a
685  positive number for fixed order, or VariableOrder. */
686  void SetOrder(int Order) const;
687 
688  virtual const FiniteElement *
690 
691  virtual int DofForGeometry(Geometry::Type GeomType) const;
692 
693  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
694  int Or) const;
695 
696  virtual const char *Name() const { return name; }
697 
698  virtual int GetContType() const { return CONTINUOUS; }
699 
701 
702  virtual ~NURBSFECollection();
703 };
704 
705 
706 /// Piecewise-(bi/tri)linear continuous finite elements.
708 {
709 private:
710  const PointFiniteElement PointFE;
711  const Linear1DFiniteElement SegmentFE;
712  const Linear2DFiniteElement TriangleFE;
713  const BiLinear2DFiniteElement QuadrilateralFE;
714  const Linear3DFiniteElement TetrahedronFE;
715  const TriLinear3DFiniteElement ParallelepipedFE;
716  const LinearWedgeFiniteElement WedgeFE;
717  const LinearPyramidFiniteElement PyramidFE;
718 public:
720 
721  virtual const FiniteElement *
723 
724  virtual int DofForGeometry(Geometry::Type GeomType) const;
725 
726  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
727  int Or) const;
728 
729  virtual const char * Name() const { return "Linear"; }
730 
731  virtual int GetContType() const { return CONTINUOUS; }
732 };
733 
734 /// Piecewise-(bi)quadratic continuous finite elements.
736 {
737 private:
738  const PointFiniteElement PointFE;
739  const Quad1DFiniteElement SegmentFE;
740  const Quad2DFiniteElement TriangleFE;
741  const BiQuad2DFiniteElement QuadrilateralFE;
742  const Quadratic3DFiniteElement TetrahedronFE;
743  const LagrangeHexFiniteElement ParallelepipedFE;
744  const H1_WedgeElement WedgeFE;
745 
746 public:
748  : FiniteElementCollection(2), ParallelepipedFE(2), WedgeFE(2) { }
749 
750  virtual const FiniteElement *
752 
753  virtual int DofForGeometry(Geometry::Type GeomType) const;
754 
755  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
756  int Or) const;
757 
758  virtual const char * Name() const { return "Quadratic"; }
759 
760  virtual int GetContType() const { return CONTINUOUS; }
761 };
762 
763 /// Version of QuadraticFECollection with positive basis functions.
765 {
766 private:
767  const QuadPos1DFiniteElement SegmentFE;
768  const BiQuadPos2DFiniteElement QuadrilateralFE;
769 
770 public:
772 
773  virtual const FiniteElement *
775 
776  virtual int DofForGeometry(Geometry::Type GeomType) const;
777 
778  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
779  int Or) const;
780 
781  virtual const char * Name() const { return "QuadraticPos"; }
782 
783  virtual int GetContType() const { return CONTINUOUS; }
784 };
785 
786 /// Piecewise-(bi)cubic continuous finite elements.
788 {
789 private:
790  const PointFiniteElement PointFE;
791  const Cubic1DFiniteElement SegmentFE;
792  const Cubic2DFiniteElement TriangleFE;
793  const BiCubic2DFiniteElement QuadrilateralFE;
794  const Cubic3DFiniteElement TetrahedronFE;
795  const LagrangeHexFiniteElement ParallelepipedFE;
796  const H1_WedgeElement WedgeFE;
797 
798 public:
801  ParallelepipedFE(3), WedgeFE(3, BasisType::ClosedUniform)
802  { }
803 
804  virtual const FiniteElement *
806 
807  virtual int DofForGeometry(Geometry::Type GeomType) const;
808 
809  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
810  int Or) const;
811 
812  virtual const char * Name() const { return "Cubic"; }
813 
814  virtual int GetContType() const { return CONTINUOUS; }
815 };
816 
817 /// Crouzeix-Raviart nonconforming elements in 2D.
819 {
820 private:
821  const P0SegmentFiniteElement SegmentFE;
822  const CrouzeixRaviartFiniteElement TriangleFE;
823  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
824 public:
826 
827  virtual const FiniteElement *
829 
830  virtual int DofForGeometry(Geometry::Type GeomType) const;
831 
832  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
833  int Or) const;
834 
835  virtual const char * Name() const { return "CrouzeixRaviart"; }
836 
837  virtual int GetContType() const { return DISCONTINUOUS; }
838 };
839 
840 /// Piecewise-linear nonconforming finite elements in 3D.
842 {
843 private:
844  const P0TriangleFiniteElement TriangleFE;
845  const P1TetNonConfFiniteElement TetrahedronFE;
846  const P0QuadFiniteElement QuadrilateralFE;
847  const RotTriLinearHexFiniteElement ParallelepipedFE;
848 
849 public:
851 
852  virtual const FiniteElement *
854 
855  virtual int DofForGeometry(Geometry::Type GeomType) const;
856 
857  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
858  int Or) const;
859 
860  virtual const char * Name() const { return "LinearNonConf3D"; }
861 
862  virtual int GetContType() const { return DISCONTINUOUS; }
863 };
864 
865 
866 /** @brief First order Raviart-Thomas finite elements in 2D. This class is kept
867  only for backward compatibility, consider using RT_FECollection instead. */
869 {
870 private:
871  const P0SegmentFiniteElement SegmentFE; // normal component on edge
872  const RT0TriangleFiniteElement TriangleFE;
873  const RT0QuadFiniteElement QuadrilateralFE;
874 public:
876 
877  virtual const FiniteElement *
879 
880  virtual int DofForGeometry(Geometry::Type GeomType) const;
881 
882  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
883  int Or) const;
884 
885  virtual const char * Name() const { return "RT0_2D"; }
886 
887  virtual int GetContType() const { return NORMAL; }
888 };
889 
890 /** @brief Second order Raviart-Thomas finite elements in 2D. This class is kept
891  only for backward compatibility, consider using RT_FECollection instead. */
893 {
894 private:
895  const P1SegmentFiniteElement SegmentFE; // normal component on edge
896  const RT1TriangleFiniteElement TriangleFE;
897  const RT1QuadFiniteElement QuadrilateralFE;
898 public:
900 
901  virtual const FiniteElement *
903 
904  virtual int DofForGeometry(Geometry::Type GeomType) const;
905 
906  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
907  int Or) const;
908 
909  virtual const char * Name() const { return "RT1_2D"; }
910 
911  virtual int GetContType() const { return NORMAL; }
912 };
913 
914 /** @brief Third order Raviart-Thomas finite elements in 2D. This class is kept
915  only for backward compatibility, consider using RT_FECollection instead. */
917 {
918 private:
919  const P2SegmentFiniteElement SegmentFE; // normal component on edge
920  const RT2TriangleFiniteElement TriangleFE;
921  const RT2QuadFiniteElement QuadrilateralFE;
922 public:
924 
925  virtual const FiniteElement *
927 
928  virtual int DofForGeometry(Geometry::Type GeomType) const;
929 
930  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
931  int Or) const;
932 
933  virtual const char * Name() const { return "RT2_2D"; }
934 
935  virtual int GetContType() const { return NORMAL; }
936 };
937 
938 /** @brief Piecewise-constant discontinuous finite elements in 2D. This class is
939  kept only for backward compatibility, consider using L2_FECollection
940  instead. */
942 {
943 private:
944  const P0TriangleFiniteElement TriangleFE;
945  const P0QuadFiniteElement QuadrilateralFE;
946 public:
948 
949  virtual const FiniteElement *
951 
952  virtual int DofForGeometry(Geometry::Type GeomType) const;
953 
954  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
955  int Or) const;
956 
957  virtual const char * Name() const { return "Const2D"; }
958 
959  virtual int GetContType() const { return DISCONTINUOUS; }
960 };
961 
962 /** @brief Piecewise-linear discontinuous finite elements in 2D. This class is
963  kept only for backward compatibility, consider using L2_FECollection
964  instead. */
966 {
967 private:
968  const Linear2DFiniteElement TriangleFE;
969  const BiLinear2DFiniteElement QuadrilateralFE;
970 
971 public:
973 
974  virtual const FiniteElement *
976 
977  virtual int DofForGeometry(Geometry::Type GeomType) const;
978 
979  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
980  int Or) const;
981 
982  virtual const char * Name() const { return "LinearDiscont2D"; }
983 
984  virtual int GetContType() const { return DISCONTINUOUS; }
985 };
986 
987 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
989 {
990 private:
991  // const CrouzeixRaviartFiniteElement TriangleFE;
992  const GaussLinear2DFiniteElement TriangleFE;
993  const GaussBiLinear2DFiniteElement QuadrilateralFE;
994 
995 public:
997 
998  virtual const FiniteElement *
1000 
1001  virtual int DofForGeometry(Geometry::Type GeomType) const;
1002 
1003  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1004  int Or) const;
1005 
1006  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
1007 
1008  virtual int GetContType() const { return DISCONTINUOUS; }
1009 };
1010 
1011 /// Linear (P1) finite elements on quadrilaterals.
1013 {
1014 private:
1015  const P1OnQuadFiniteElement QuadrilateralFE;
1016 public:
1018  virtual const FiniteElement *
1019  FiniteElementForGeometry(Geometry::Type GeomType) const;
1020  virtual int DofForGeometry(Geometry::Type GeomType) const;
1021  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1022  int Or) const;
1023  virtual const char * Name() const { return "P1OnQuad"; }
1024  virtual int GetContType() const { return DISCONTINUOUS; }
1025 };
1026 
1027 /** @brief Piecewise-quadratic discontinuous finite elements in 2D. This class
1028  is kept only for backward compatibility, consider using L2_FECollection
1029  instead. */
1031 {
1032 private:
1033  const Quad2DFiniteElement TriangleFE;
1034  const BiQuad2DFiniteElement QuadrilateralFE;
1035 
1036 public:
1038 
1039  virtual const FiniteElement *
1040  FiniteElementForGeometry(Geometry::Type GeomType) const;
1041 
1042  virtual int DofForGeometry(Geometry::Type GeomType) const;
1043 
1044  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1045  int Or) const;
1046 
1047  virtual const char * Name() const { return "QuadraticDiscont2D"; }
1048  virtual int GetContType() const { return DISCONTINUOUS; }
1049 };
1050 
1051 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
1053 {
1054 private:
1055  const BiQuadPos2DFiniteElement QuadrilateralFE;
1056 
1057 public:
1059  virtual const FiniteElement *
1060  FiniteElementForGeometry(Geometry::Type GeomType) const;
1061  virtual int DofForGeometry(Geometry::Type GeomType) const;
1062  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1063  int Or) const
1064  { return NULL; }
1065  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
1066  virtual int GetContType() const { return DISCONTINUOUS; }
1067 };
1068 
1069 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
1071 {
1072 private:
1073  // const Quad2DFiniteElement TriangleFE;
1074  const GaussQuad2DFiniteElement TriangleFE;
1075  const GaussBiQuad2DFiniteElement QuadrilateralFE;
1076 
1077 public:
1079 
1080  virtual const FiniteElement *
1081  FiniteElementForGeometry(Geometry::Type GeomType) const;
1082 
1083  virtual int DofForGeometry(Geometry::Type GeomType) const;
1084 
1085  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1086  int Or) const;
1087 
1088  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
1089  virtual int GetContType() const { return DISCONTINUOUS; }
1090 };
1091 
1092 /** @brief Piecewise-cubic discontinuous finite elements in 2D. This class is
1093  kept only for backward compatibility, consider using L2_FECollection
1094  instead. */
1096 {
1097 private:
1098  const Cubic2DFiniteElement TriangleFE;
1099  const BiCubic2DFiniteElement QuadrilateralFE;
1100 
1101 public:
1103 
1104  virtual const FiniteElement *
1105  FiniteElementForGeometry(Geometry::Type GeomType) const;
1106 
1107  virtual int DofForGeometry(Geometry::Type GeomType) const;
1108 
1109  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1110  int Or) const;
1111 
1112  virtual const char * Name() const { return "CubicDiscont2D"; }
1113  virtual int GetContType() const { return DISCONTINUOUS; }
1114 };
1115 
1116 /** @brief Piecewise-constant discontinuous finite elements in 3D. This class is
1117  kept only for backward compatibility, consider using L2_FECollection
1118  instead. */
1120 {
1121 private:
1122  const P0TetFiniteElement TetrahedronFE;
1123  const P0HexFiniteElement ParallelepipedFE;
1124  const P0WdgFiniteElement WedgeFE;
1125  const P0PyrFiniteElement PyramidFE;
1126 
1127 public:
1129 
1130  virtual const FiniteElement *
1131  FiniteElementForGeometry(Geometry::Type GeomType) const;
1132 
1133  virtual int DofForGeometry(Geometry::Type GeomType) const;
1134 
1135  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1136  int Or) const;
1137 
1138  virtual const char * Name() const { return "Const3D"; }
1139  virtual int GetContType() const { return DISCONTINUOUS; }
1140 };
1141 
1142 /** @brief Piecewise-linear discontinuous finite elements in 3D. This class is
1143  kept only for backward compatibility, consider using L2_FECollection
1144  instead. */
1146 {
1147 private:
1148  const Linear3DFiniteElement TetrahedronFE;
1149  const LinearPyramidFiniteElement PyramidFE;
1150  const LinearWedgeFiniteElement WedgeFE;
1151  const TriLinear3DFiniteElement ParallelepipedFE;
1152 
1153 public:
1155 
1156  virtual const FiniteElement *
1157  FiniteElementForGeometry(Geometry::Type GeomType) const;
1158 
1159  virtual int DofForGeometry(Geometry::Type GeomType) const;
1160 
1161  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1162  int Or) const;
1163 
1164  virtual const char * Name() const { return "LinearDiscont3D"; }
1165  virtual int GetContType() const { return DISCONTINUOUS; }
1166 };
1167 
1168 /** @brief Piecewise-quadratic discontinuous finite elements in 3D. This class
1169  is kept only for backward compatibility, consider using L2_FECollection
1170  instead. */
1172 {
1173 private:
1174  const Quadratic3DFiniteElement TetrahedronFE;
1175  const LagrangeHexFiniteElement ParallelepipedFE;
1176 
1177 public:
1179  : FiniteElementCollection(2), ParallelepipedFE(2) { }
1180 
1181  virtual const FiniteElement *
1182  FiniteElementForGeometry(Geometry::Type GeomType) const;
1183 
1184  virtual int DofForGeometry(Geometry::Type GeomType) const;
1185 
1186  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1187  int Or) const;
1188 
1189  virtual const char * Name() const { return "QuadraticDiscont3D"; }
1190  virtual int GetContType() const { return DISCONTINUOUS; }
1191 };
1192 
1193 /// Finite element collection on a macro-element.
1195 {
1196 private:
1197  const PointFiniteElement PointFE;
1198  const RefinedLinear1DFiniteElement SegmentFE;
1199  const RefinedLinear2DFiniteElement TriangleFE;
1200  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
1201  const RefinedLinear3DFiniteElement TetrahedronFE;
1202  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
1203 
1204 public:
1206 
1207  virtual const FiniteElement *
1208  FiniteElementForGeometry(Geometry::Type GeomType) const;
1209 
1210  virtual int DofForGeometry(Geometry::Type GeomType) const;
1211 
1212  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1213  int Or) const;
1214 
1215  virtual const char * Name() const { return "RefinedLinear"; }
1216  virtual int GetContType() const { return CONTINUOUS; }
1217 };
1218 
1219 /** @brief Lowest order Nedelec finite elements in 3D. This class is kept only
1220  for backward compatibility, consider using the new ND_FECollection
1221  instead. */
1223 {
1224 private:
1225  const Nedelec1HexFiniteElement HexahedronFE;
1226  const Nedelec1TetFiniteElement TetrahedronFE;
1227  const Nedelec1WdgFiniteElement WedgeFE;
1228  const Nedelec1PyrFiniteElement PyramidFE;
1229 
1230 public:
1232 
1233  virtual const FiniteElement *
1234  FiniteElementForGeometry(Geometry::Type GeomType) const;
1235 
1236  virtual int DofForGeometry(Geometry::Type GeomType) const;
1237 
1238  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1239  int Or) const;
1240 
1241  virtual const char * Name() const { return "ND1_3D"; }
1242  virtual int GetContType() const { return TANGENTIAL; }
1243 };
1244 
1245 /** @brief First order Raviart-Thomas finite elements in 3D. This class is kept
1246  only for backward compatibility, consider using RT_FECollection instead. */
1248 {
1249 private:
1250  const P0TriangleFiniteElement TriangleFE;
1251  const P0QuadFiniteElement QuadrilateralFE;
1252  const RT0HexFiniteElement HexahedronFE;
1253  const RT0TetFiniteElement TetrahedronFE;
1254  const RT0WdgFiniteElement WedgeFE;
1255  const RT0PyrFiniteElement PyramidFE;
1256 public:
1258 
1259  virtual const FiniteElement *
1260  FiniteElementForGeometry(Geometry::Type GeomType) const;
1261 
1262  virtual int DofForGeometry(Geometry::Type GeomType) const;
1263 
1264  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1265  int Or) const;
1266 
1267  virtual const char * Name() const { return "RT0_3D"; }
1268  virtual int GetContType() const { return NORMAL; }
1269 };
1270 
1271 /** @brief Second order Raviart-Thomas finite elements in 3D. This class is kept
1272  only for backward compatibility, consider using RT_FECollection instead. */
1274 {
1275 private:
1276  const Linear2DFiniteElement TriangleFE;
1277  const BiLinear2DFiniteElement QuadrilateralFE;
1278  const RT1HexFiniteElement HexahedronFE;
1279 public:
1281 
1282  virtual const FiniteElement *
1283  FiniteElementForGeometry(Geometry::Type GeomType) const;
1284 
1285  virtual int DofForGeometry(Geometry::Type GeomType) const;
1286 
1287  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1288  int Or) const;
1289 
1290  virtual const char * Name() const { return "RT1_3D"; }
1291  virtual int GetContType() const { return NORMAL; }
1292 };
1293 
1294 /// Discontinuous collection defined locally by a given finite element.
1296 {
1297 private:
1298  char d_name[32];
1299  Geometry::Type GeomType;
1300  FiniteElement *Local_Element;
1301 
1302 public:
1303  Local_FECollection(const char *fe_name);
1304 
1306  Geometry::Type GeomType_) const
1307  { return (GeomType == GeomType_) ? Local_Element : NULL; }
1308  virtual int DofForGeometry(Geometry::Type GeomType_) const
1309  { return (GeomType == GeomType_) ? Local_Element->GetDof() : 0; }
1310  virtual const int *DofOrderForOrientation(Geometry::Type GeomType_,
1311  int Or) const
1312  { return NULL; }
1313  virtual const char *Name() const { return d_name; }
1314 
1315  virtual ~Local_FECollection() { delete Local_Element; }
1316  virtual int GetContType() const { return DISCONTINUOUS; }
1317 };
1318 
1319 }
1320 
1321 #endif
Abstract class for all finite elements.
Definition: fe_base.hpp:233
virtual const char * Name() const
Definition: fe_coll.hpp:581
A 2D 2nd order Raviart-Thomas vector element on a triangle.
virtual const char * Name() const
Definition: fe_coll.hpp:418
virtual int GetContType() const
Definition: fe_coll.hpp:631
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:649
virtual const char * Name() const
Definition: fe_coll.hpp:526
virtual int GetContType() const
Definition: fe_coll.hpp:911
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:2045
virtual int GetContType() const
Definition: fe_coll.hpp:959
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1176
Piecewise-cubic discontinuous finite elements in 2D. This class is kept only for backward compatibili...
Definition: fe_coll.hpp:1095
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:1070
virtual const char * Name() const
Definition: fe_coll.hpp:860
BiLinear2DFiniteElement QuadrilateralFE
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:388
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:461
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:462
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:414
A 2D 1st order Raviart-Thomas vector element on a square.
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3095
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1462
A 1D constant element on a segment.
A 1D linear element with nodes at 1/3 and 2/3 (trace of RT1)
A 1D quadratic element with nodes at the Gaussian points (trace of RT2)
virtual const char * Name() const
Definition: fe_coll.hpp:1267
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1420
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:539
virtual const char * Name() const
Definition: fe_coll.hpp:1047
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1084
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:3483
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:870
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:1618
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:2453
A 3D 1st order Nedelec element on a tetrahedron.
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:983
virtual const char * Name() const
Definition: fe_coll.hpp:933
ErrorMode
How to treat errors in FiniteElementForGeometry() calls.
Definition: fe_coll.hpp:245
Linear1DFiniteElement SegmentFE
Definition: segment.cpp:49
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1603
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:3193
Field is discontinuous across element interfaces.
Definition: fe_coll.hpp:48
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:884
static const int NumGeom
Definition: geom.hpp:42
A 3D 1st order Nedelec element on a pyramid.
virtual int GetContType() const =0
const FiniteElement * GetFE(Geometry::Type geom, int p) const
Variable order version of FiniteElementForGeometry().
Definition: fe_coll.hpp:191
virtual int GetContType() const
Definition: fe_coll.hpp:837
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:927
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2874
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:1437
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1488
A 2D 3rd order Raviart-Thomas vector element on a square.
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
Definition: fe_coll.hpp:507
virtual ~L2_FECollection()
Definition: fe_coll.cpp:2366
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:1478
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:1356
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:679
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition: fe_coll.hpp:225
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const =0
int GetClosedBasisType() const
Definition: fe_coll.hpp:422
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1951
virtual const char * Name() const
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:3360
Piecewise-(bi/tri)linear continuous finite elements.
Definition: fe_coll.hpp:707
virtual int GetContType() const
Definition: fe_coll.hpp:527
virtual int GetContType() const
Definition: fe_coll.hpp:419
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:2603
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:3020
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:634
A 2D quadratic element on triangle with nodes at the vertices and midpoints of the triangle...
virtual const char * Name() const
Definition: fe_coll.hpp:1112
A 0D point finite element.
Arbitrary order 3D H(curl)-trace finite elements in 2D defined on the interface between mesh elements...
Definition: fe_coll.hpp:591
A 3D 0th order Raviert-Thomas element on a cube.
virtual const char * Name() const
Definition: fe_coll.hpp:812
A 1D quadratic finite element with uniformly spaced nodes.
virtual int GetContType() const
Definition: fe_coll.hpp:1242
virtual int GetContType() const
Definition: fe_coll.hpp:760
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
Definition: fe_coll.hpp:561
int GetClosedBasisType() const
Definition: fe_coll.hpp:486
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:1237
virtual ~RT_R2D_FECollection()
Definition: fe_coll.cpp:3386
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1069
int GetBasisType() const
Definition: fe_coll.hpp:282
const int base_p
Order as returned by GetOrder().
Definition: fe_coll.hpp:235
Second order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:892
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:2468
virtual int GetContType() const
Definition: fe_coll.hpp:935
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:362
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:1062
A 3D constant element on a tetrahedron.
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:522
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:2670
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1553
RT_R2D_FECollection(const int p, const int dim, const int map_type, const bool signs, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3307
A 3D Crouzeix-Raviart element on the tetrahedron.
virtual const char * Name() const
Definition: fe_coll.hpp:1138
virtual const char * Name() const
Definition: fe_coll.hpp:758
A 2D bi-cubic element on a square with uniformly spaces nodes.
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
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:370
virtual const char * Name() const
Definition: fe_coll.hpp:957
Finite element collection on a macro-element.
Definition: fe_coll.hpp:1194
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:1100
A 2D refined bi-linear FE on a square.
Possible basis types. Note that not all elements can use all BasisType(s).
Definition: fe_base.hpp:25
virtual const char * Name() const
Definition: fe_coll.hpp:1023
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1118
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:1313
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default)...
Definition: fe_coll.cpp:3458
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1406
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:779
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:738
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:473
std::function< double(const Vector &)> f(double mass_coeff)
Definition: lor_mms.hpp:30
virtual int GetRangeType(int dim) const
Definition: fe_coll.cpp:40
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2638
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:2065
virtual const FiniteElement * FiniteElementForDim(int dim) const
Returns the first non-NULL FiniteElement for the given dimension.
Definition: fe_coll.cpp:26
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.cpp:361
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:787
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:3089
Class for finite elements with basis functions that return scalar values.
Definition: fe_base.hpp:649
virtual int GetContType() const
Definition: fe_coll.hpp:1165
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:549
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1988
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2332
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1323
A 2D bi-linear element on a square with nodes at the "Gaussian" points.
A 1D refined linear element.
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:266
A 2D refined linear element on a triangle.
A 2D refined linear element on a tetrahedron.
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:417
const int * GetDofOrdering(Geometry::Type geom, int p, int ori) const
Variable order version of DofOrderForOrientation().
Definition: fe_coll.hpp:215
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:387
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:1310
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1297
virtual int GetContType() const
Definition: fe_coll.hpp:483
A 2D constant element on a triangle.
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:387
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:291
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:511
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:969
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1205
virtual int GetDerivType(int dim) const
Definition: fe_coll.cpp:70
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:860
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:913
A linear element defined on a triangular prism.
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:653
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1221
Piecewise-linear discontinuous finite elements in 2D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:965
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:2902
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1139
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:565
virtual const char * Name() const
Definition: fe_coll.hpp:553
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1640
virtual const char * Name() const
Definition: fe_coll.hpp:630
virtual int GetContType() const
Definition: fe_coll.hpp:1268
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:2952
A 3D 0th order Raviert-Thomas element on a wedge.
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:319
virtual int GetContType() const
Definition: fe_coll.hpp:360
ND_R2D_Trace_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3233
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2690
ErrorMode error_mode
How to treat errors in FiniteElementForGeometry() calls.
Definition: fe_coll.hpp:255
A 3D refined tri-linear element on a cube.
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:1295
Arbitrary order H(curl)-trace finite elements defined on the interface between mesh elements (faces...
Definition: fe_coll.hpp:498
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:1052
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2940
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:941
A 3D quadratic element on a tetrahedron with uniformly spaced nodes.
RT_R1D_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3040
virtual int GetContType() const
Definition: fe_coll.hpp:698
virtual int GetContType() const
Definition: fe_coll.hpp:1216
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1536
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:110
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1108
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:697
Array< FiniteElementCollection * > var_orders
Definition: fe_coll.hpp:242
Third order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:916
First order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:1247
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:623
virtual const char * Name() const
Definition: fe_coll.hpp:80
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:275
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:3492
virtual int GetContType() const
Definition: fe_coll.hpp:984
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:764
virtual ~Local_FECollection()
Definition: fe_coll.hpp:1315
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1338
virtual ~H1_FECollection()
Definition: fe_coll.cpp:2032
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:1119
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:1025
virtual const char * Name() const
Definition: fe_coll.hpp:909
A 3D linear element on a tetrahedron with nodes at the vertices of the tetrahedron.
A 1D linear element with nodes on the endpoints.
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:898
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:1521
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:512
A 2D Crouzeix-Raviart finite element on square.
A 2D linear element on triangle with nodes at the vertices of the triangle.
void Reset() const
Definition: fe_coll.hpp:672
A quadratic element on triangle with nodes at the "Gaussian" points.
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:941
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:265
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:841
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1189
virtual int GetContType() const
Definition: fe_coll.hpp:814
Tangential components of vector field.
Definition: fe_coll.hpp:46
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:1273
virtual int GetContType() const
Definition: fe_coll.hpp:731
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:818
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:605
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1503
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:380
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:377
virtual const char * Name() const
Definition: fe_coll.hpp:1313
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:302
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:3513
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:397
H1Ser_FECollection(const int p, const int dim=2)
Definition: fe_coll.hpp:312
A 2D 2nd order Raviart-Thomas vector element on a square.
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]...
A 3D 1st order Raviert-Thomas element on a cube.
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1380
A 2D bi-quadratic element on a square with uniformly spaced nodes.
int GetNumDof(Geometry::Type geom, int p) const
Variable order version of DofForGeometry().
Definition: fe_coll.hpp:203
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2588
virtual const char * Name() const
Definition: fe_coll.hpp:280
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:425
A 3D tri-linear element on a cube with nodes at the vertices of the cube.
ND_R1D_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2971
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:955
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3520
A 3D constant element on a wedge.
virtual const char * Name() const
Definition: fe_coll.hpp:1065
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:97
virtual int GetContType() const
Definition: fe_coll.hpp:1024
Piecewise-quadratic discontinuous finite elements in 3D. This class is kept only for backward compati...
Definition: fe_coll.hpp:1171
virtual int GetContType() const
Definition: fe_coll.hpp:281
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3370
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
virtual int GetDerivMapType(int dim) const
Definition: fe_coll.cpp:80
A 3D constant element on a cube.
virtual int GetContType() const
Definition: fe_coll.hpp:1139
A 2D bi-linear element on a square with nodes at the vertices of the square.
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:988
ND_R2D_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3110
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:116
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1033
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType_) const
Definition: fe_coll.hpp:1305
A linear element on a triangle with nodes at the 3 "Gaussian" points.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:831
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:1222
Arbitrary order H1-conforming (continuous) finite elements with positive basis functions.
Definition: fe_coll.hpp:299
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3203
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1589
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1153
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2920
int HasFaceDofs(Geometry::Type geom, int p) const
Definition: fe_coll.cpp:90
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:682
A 2D bi-quadratic element on a square with nodes at the 9 "Gaussian" points.
virtual const char * Name() const
Definition: fe_coll.hpp:1215
virtual int GetMapType(int dim) const
Definition: fe_coll.cpp:60
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:1131
MFEM_EXPORT class LinearWedgeFiniteElement WedgeFE
Definition: fe.cpp:40
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:751
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1364
virtual const char * Name() const
Definition: fe_coll.hpp:885
virtual int GetContType() const
Definition: fe_coll.hpp:582
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:1061
A 1D quadratic positive element utilizing the 2nd order Bernstein basis.
Definition: fe_pos.hpp:107
A 3D 1st order Nedelec element on a wedge.
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:714
int dim
Definition: ex24.cpp:53
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:3321
Piecewise-linear discontinuous finite elements in 3D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:1145
RT_R2D_Trace_FECollection(const int p, const int dim, const int map_type=FiniteElement::INTEGRAL, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:3396
virtual int DofForGeometry(Geometry::Type GeomType_) const
Definition: fe_coll.hpp:1308
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:724
virtual int GetContType() const
Definition: fe_coll.hpp:783
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:566
virtual const char * Name() const
Definition: fe_coll.hpp:358
virtual int GetContType() const
Definition: fe_coll.hpp:554
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:606
virtual const char * Name() const
Definition: fe_coll.hpp:1006
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:1398
int GetOpenBasisType() const
Definition: fe_coll.hpp:423
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:2650
Tensor products of 1D Lagrange1DFiniteElement (only degree 2 is functional)
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:574
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:626
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:735
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1444
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:577
virtual const char * Name() const
Definition: fe_coll.hpp:982
First order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:868
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:1966
A 3D constant element on a pyramid.
A 2D constant element on a square.
virtual ~RT_R1D_FECollection()
Definition: fe_coll.cpp:3101
A linear element defined on a square pyramid.
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:2007
virtual int GetDerivRangeType(int dim) const
Definition: fe_coll.cpp:50
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:846
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:454
virtual const char * Name() const
Definition: fe_coll.hpp:1088
virtual const char * Name() const
Definition: fe_coll.hpp:835
A 2D cubic element on a triangle with uniformly spaced nodes.
A 1D cubic element with uniformly spaced nodes.
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:434
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:3507
A 2D 1st order Raviart-Thomas vector element on a triangle.
Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 2D.
Definition: fe_coll.hpp:600
virtual ~ND_R1D_FECollection()
Definition: fe_coll.cpp:3031
virtual int GetContType() const
Definition: fe_coll.hpp:1113
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
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:2347
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:671
Second order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:1273
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1011
A 2D linear element on a square with 3 nodes at the vertices of the lower left triangle.
virtual StatelessDofTransformation * DofTransformationForGeometry(Geometry::Type GeomType) const
Returns a DoF transformation object compatible with this basis and geometry type. ...
Definition: fe_coll.hpp:68
A 3D 0th order Raviert-Thomas element on a pyramid.
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3026
Field is continuous across element interfaces.
Definition: fe_coll.hpp:45
virtual const char * Name() const
Definition: fe_coll.hpp:482
virtual int GetContType() const
Definition: fe_coll.hpp:862
Piecewise-quadratic discontinuous finite elements in 2D. This class is kept only for backward compati...
Definition: fe_coll.hpp:1030
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1047
Arbitrary order H1 elements in 3D on a wedge.
Definition: fe_h1.hpp:130
A 3D 0th order Raviert-Thomas element on a tetrahedron.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1245
virtual const char * Name() const
Definition: fe_coll.hpp:1189
int GetBasisType() const
Definition: fe_coll.hpp:368
Arbitrary order 3D "H^{-1/2}-conforming" face finite elements defined on the interface between mesh e...
Definition: fe_coll.hpp:640
virtual const char * Name() const
Definition: fe_coll.hpp:1164
A 2D 3rd order Raviart-Thomas vector element on a triangle.
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:3416
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2621
virtual int GetContType() const
Definition: fe_coll.hpp:1316
virtual int GetContType() const
Definition: fe_coll.hpp:1291
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1281
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:3470
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:1012
void InitVarOrder(int p) const
Definition: fe_coll.cpp:369
Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 1D.
Definition: fe_coll.hpp:534
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:761
virtual const char * Name() const
Definition: fe_coll.hpp:781
virtual int GetContType() const
Definition: fe_coll.hpp:887
An arbitrary order 3D NURBS element on a cube.
Definition: fe_nurbs.hpp:119
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:489
An arbitrary order 1D NURBS element on a segment.
Definition: fe_nurbs.hpp:67
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:546
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:462
virtual const char * Name() const
Definition: fe_coll.hpp:1241
virtual const char * Name() const
Definition: fe_coll.hpp:696
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:538
An arbitrary order 2D NURBS element on a square.
Definition: fe_nurbs.hpp:87
virtual StatelessDofTransformation * DofTransformationForGeometry(Geometry::Type GeomType) const
Returns a DoF transformation object compatible with this basis and geometry type. ...
Definition: fe_coll.cpp:2890
A 3D 1st order Nedelec element on a cube.
virtual const char * Name() const
Definition: fe_coll.hpp:1290
A 2D Crouzeix-Raviart element on triangle.
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:1571
virtual ~ND_R2D_FECollection()
Definition: fe_coll.cpp:3223
int GetOpenBasisType() const
Definition: fe_coll.hpp:487
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:1168
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:327
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:348
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:998
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:519
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1259