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