MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
fe_coll.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_FE_COLLECTION
13 #define MFEM_FE_COLLECTION
14 
15 #include "../config/config.hpp"
16 #include "geom.hpp"
17 #include "fe.hpp"
18 
19 namespace mfem
20 {
21 
22 /** @brief Collection of finite elements from the same family in multiple
23  dimensions. This class is used to match the degrees of freedom of a
24  FiniteElementSpace between elements, and to provide the finite element
25  restriction from an element to its boundary. */
27 {
28 protected:
29  template <Geometry::Type geom>
30  static inline void GetNVE(int &nv, int &ne);
31 
32  template <Geometry::Type geom, typename v_t>
33  static inline void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo,
34  const int edge_info);
35 
36  template <Geometry::Type geom, Geometry::Type f_geom,
37  typename v_t, typename e_t, typename eo_t>
38  static inline void GetFace(int &nv, v_t &v, int &ne, e_t &e, eo_t &eo,
39  int &nf, int &f, Geometry::Type &fg, int &fo,
40  const int face_info);
41 
42 public:
43  /** @brief Enumeration for ContType: defines the continuity of the field
44  across element interfaces. */
45  enum { CONTINUOUS, ///< Field is continuous across element interfaces
46  TANGENTIAL, ///< Tangential components of vector field
47  NORMAL, ///< Normal component of vector field
48  DISCONTINUOUS ///< Field is discontinuous across element interfaces
49  };
50 
51  virtual const FiniteElement *
52  FiniteElementForGeometry(Geometry::Type GeomType) const = 0;
53 
54  virtual int DofForGeometry(Geometry::Type GeomType) const = 0;
55 
56  /** @brief Returns an array, say p, that maps a local permuted index i to a
57  local base index: base_i = p[i].
58 
59  @note Only provides information about interior dofs. See
60  FiniteElementCollection::SubDofOrder if interior \a and boundary dof
61  order is needed. */
62  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
63  int Or) const = 0;
64 
65  virtual const char * Name() const { return "Undefined"; }
66 
67  virtual int GetContType() const = 0;
68 
69  int HasFaceDofs(Geometry::Type geom, int p) const;
70 
72  Geometry::Type GeomType) const
73  {
74  return FiniteElementForGeometry(GeomType);
75  }
76 
78 
79  virtual ~FiniteElementCollection();
80 
81  /** @brief Factory method: return a newly allocated FiniteElementCollection
82  according to the given name. */
83  /**
84  | FEC Name | Space | Order | BasisType | FiniteElement::MapT | Notes |
85  | :------: | :---: | :---: | :-------: | :-----: | :---: |
86  | H1_[DIM]_[ORDER] | H1 | * | 1 | VALUE | H1 nodal elements |
87  | H1@[BTYPE]_[DIM]_[ORDER] | H1 | * | * | VALUE | H1 nodal elements |
88  | H1Pos_[DIM]_[ORDER] | H1 | * | 1 | VALUE | H1 nodal elements |
89  | 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) |
90  | 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) |
91  | 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) |
92  | ND_[DIM]_[ORDER] | H(curl) | * | 1 / 0 | H_CURL | Nedelec vector elements |
93  | ND@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(curl) | * | * / * | H_CURL | Nedelec vector elements |
94  | 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) |
95  | 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) |
96  | RT_[DIM]_[ORDER] | H(div) | * | 1 / 0 | H_DIV | Raviart-Thomas vector elements |
97  | RT@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(div) | * | * / * | H_DIV | Raviart-Thomas vector elements |
98  | 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) |
99  | 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) |
100  | 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) |
101  | 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) |
102  | L2_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinuous L2 elements |
103  | L2_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinuous L2 elements |
104  | L2Int_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinuous L2 elements |
105  | L2Int_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinuous L2 elements |
106  | DG_Iface_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
107  | DG_Iface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
108  | DG_IntIface_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
109  | DG_IntIface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
110  | NURBS[ORDER] | - | * | - | VALUE | Non-Uniform Rational B-Splines (NURBS) elements |
111  | LinearNonConf3D | - | 1 | 1 | VALUE | Piecewise-linear nonconforming finite elements in 3D |
112  | CrouzeixRaviart | - | - | - | - | Crouzeix-Raviart nonconforming elements in 2D |
113  | Local_[FENAME] | - | - | - | - | Special collection that builds a local version out of the FENAME collection |
114  |-|-|-|-|-|-|
115  | Linear | H1 | 1 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
116  | Quadratic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
117  | QuadraticPos | H1 | 2 | 2 | VALUE | Left in for backward compatibility, consider using H1_ |
118  | Cubic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
119  | Const2D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
120  | Const3D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
121  | LinearDiscont2D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
122  | GaussLinearDiscont2D | L2 | 1 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
123  | P1OnQuad | H1 | 1 | 1 | VALUE | Linear P1 element with 3 nodes on a square |
124  | QuadraticDiscont2D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
125  | QuadraticPosDiscont2D | L2 | 2 | 2 | VALUE | Left in for backward compatibility, consider using L2_ |
126  | GaussQuadraticDiscont2D | L2 | 2 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
127  | CubicDiscont2D | L2 | 3 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
128  | LinearDiscont3D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
129  | QuadraticDiscont3D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
130  | ND1_3D | H(Curl) | 1 | 1 / 0 | H_CURL | Left in for backward compatibility, consider using ND_ |
131  | RT0_2D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
132  | RT1_2D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
133  | RT2_2D | H(Div) | 3 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
134  | RT0_3D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
135  | RT1_3D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
136 
137  | Tag | Description |
138  | :------: | :--------: |
139  | [DIM] | Dimension of the elements (1D, 2D, 3D) |
140  | [ORDER] | Approximation order of the elements (P0, P1, P2, ...) |
141  | [BTYPE] | BasisType of the element (0-GaussLegendre, 1 - GaussLobatto, 2-Bernstein, 3-OpenUniform, 4-CloseUniform, 5-OpenHalfUniform) |
142  | [OBTYPE] | Open BasisType of the element for elements which have both types |
143  | [CBTYPE] | Closed BasisType of the element for elements which have both types |
144 
145  [FENAME] Is a special case for the Local FEC which generates a local version of a given
146  FEC. It is selected from one of (BiCubic2DFiniteElement, Quad_Q3, Nedelec1HexFiniteElement,
147  Hex_ND1, H1_[DIM]_[ORDER],H1Pos_[DIM]_[ORDER], L2_[DIM]_[ORDER] )
148  */
149  static FiniteElementCollection *New(const char *name);
150 
151  /** @brief Get the local dofs for a given sub-manifold.
152 
153  Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex, 1D
154  - edge, 2D - face) including those on its boundary. The local index of the
155  sub-manifold (inside Geom) and its orientation are given by the parameter
156  Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed that 0 <=
157  SDim <= Dim(Geom). */
158  void SubDofOrder(Geometry::Type Geom, int SDim, int Info,
159  Array<int> &dofs) const;
160 
161  /// Variable order version of FiniteElementForGeometry().
162  /** The order parameter @a p represents the order of the highest-dimensional
163  FiniteElement%s the fixed-order collection we want to query. In general,
164  this order is different from the order of the returned FiniteElement. */
165  const FiniteElement *GetFE(Geometry::Type geom, int p) const
166  {
167  if (p == base_p) { return FiniteElementForGeometry(geom); }
168  if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
169  return var_orders[p]->FiniteElementForGeometry(geom);
170  }
171 
172  /// Variable order version of DofForGeometry().
173  /** The order parameter @a p represents the order of the highest-dimensional
174  FiniteElement%s the fixed-order collection we want to query. In general,
175  this order is different from the order of the element corresponding to
176  @a geom in that fixed-order collection. */
177  int GetNumDof(Geometry::Type geom, int p) const
178  {
179  if (p == base_p) { return DofForGeometry(geom); }
180  if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
181  return var_orders[p]->DofForGeometry(geom);
182  }
183 
184  /// Variable order version of DofOrderForOrientation().
185  /** The order parameter @a p represents the order of the highest-dimensional
186  FiniteElement%s the fixed-order collection we want to query. In general,
187  this order is different from the order of the element corresponding to
188  @a geom in that fixed-order collection. */
189  const int *GetDofOrdering(Geometry::Type geom, int p, int ori) const
190  {
191  if (p == base_p) { return DofOrderForOrientation(geom, ori); }
192  if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
193  return var_orders[p]->DofOrderForOrientation(geom, ori);
194  }
195 
196  /** @brief Return the order (polynomial degree) of the FE collection,
197  corresponding to the order/degree returned by FiniteElement::GetOrder()
198  of the highest-dimensional FiniteElement%s defined by the collection. */
199  int GetOrder() const { return base_p; }
200 
201  /// Instantiate a new collection of the same type with a different order.
202  /** Generally, the order parameter @a p is NOT the same as the parameter @a p
203  used by some of the constructors of derived classes. Instead, this @a p
204  represents the order of the new FE collection as it will be returned by
205  its GetOrder() method. */
206  virtual FiniteElementCollection *Clone(int p) const;
207 
208 protected:
209  const int base_p; ///< Order as returned by GetOrder().
210 
213 
214  void InitVarOrder(int p) const;
215 
217 };
218 
219 /// Arbitrary order H1-conforming (continuous) finite elements.
221 {
222 
223 protected:
224  int dim, b_type;
225  char h1_name[32];
228  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8], *TetDofOrd[24];
229 
230 public:
231  explicit H1_FECollection(const int p, const int dim = 3,
232  const int btype = BasisType::GaussLobatto);
233 
235  Geometry::Type GeomType) const;
236  virtual int DofForGeometry(Geometry::Type GeomType) const
237  { return H1_dof[GeomType]; }
238  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
239  int Or) const;
240 
241  virtual const char *Name() const { return h1_name; }
242  virtual int GetContType() const { return CONTINUOUS; }
243  int GetBasisType() const { return b_type; }
244 
246 
247  /// Get the Cartesian to local H1 dof map
248  const int *GetDofMap(Geometry::Type GeomType) const;
249  /// Variable order version of GetDofMap
250  const int *GetDofMap(Geometry::Type GeomType, int p) const;
251 
253  { return new H1_FECollection(p, dim, b_type); }
254 
255  virtual ~H1_FECollection();
256 };
257 
258 /** @brief Arbitrary order H1-conforming (continuous) finite elements with
259  positive basis functions. */
261 {
262 public:
263  explicit H1Pos_FECollection(const int p, const int dim = 3)
264  : H1_FECollection(p, dim, BasisType::Positive) { }
265 };
266 
267 
268 /** Arbitrary order H1-conforming (continuous) serendipity finite elements;
269  Current implementation works in 2D only; 3D version is in development. */
271 {
272 public:
273  explicit H1Ser_FECollection(const int p, const int dim = 2)
274  : H1_FECollection(p, dim, BasisType::Serendipity) { };
275 };
276 
277 /** @brief Arbitrary order "H^{1/2}-conforming" trace finite elements defined on
278  the interface between mesh elements (faces,edges,vertices); these are the
279  trace FEs of the H1-conforming FEs. */
281 {
282 public:
283  H1_Trace_FECollection(const int p, const int dim,
284  const int btype = BasisType::GaussLobatto);
285 };
286 
287 /// Arbitrary order "L2-conforming" discontinuous finite elements.
289 {
290 private:
291  int dim;
292  int b_type; // BasisType
293  int m_type; // map type
294  char d_name[32];
297  int *SegDofOrd[2]; // for rotating segment dofs in 1D
298  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
299  int *TetDofOrd[24]; // for rotating tetrahedron dofs in 3D
300  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
301 
302 public:
303  L2_FECollection(const int p, const int dim,
304  const int btype = BasisType::GaussLegendre,
305  const int map_type = FiniteElement::VALUE);
306 
308  Geometry::Type GeomType) const;
309  virtual int DofForGeometry(Geometry::Type GeomType) const
310  {
311  if (L2_Elements[GeomType])
312  {
313  return L2_Elements[GeomType]->GetDof();
314  }
315  return 0;
316  }
317  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
318  int Or) const;
319  virtual const char *Name() const { return d_name; }
320 
321  virtual int GetContType() const { return DISCONTINUOUS; }
322 
324  Geometry::Type GeomType) const
325  {
326  return Tr_Elements[GeomType];
327  }
328 
329  int GetBasisType() const { return b_type; }
330 
332  { return new L2_FECollection(p, dim, b_type, m_type); }
333 
334  virtual ~L2_FECollection();
335 };
336 
337 /// Declare an alternative name for L2_FECollection = DG_FECollection
339 
340 /// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
342 {
343 protected:
344  int dim;
345  int cb_type; // closed BasisType
346  int ob_type; // open BasisType
347  char rt_name[32];
350  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
351 
352  // Initialize only the face elements
353  void InitFaces(const int p, const int dim, const int map_type,
354  const bool signs);
355 
356  // Constructor used by the constructor of the RT_Trace_FECollection and
357  // DG_Interface_FECollection classes
358  RT_FECollection(const int p, const int dim, const int map_type,
359  const bool signs,
360  const int ob_type = BasisType::GaussLegendre);
361 
362 public:
363  /// Construct an H(div)-conforming Raviart-Thomas FE collection, RT_p.
364  /** The index @a p corresponds to the space RT_p, as typically denoted in the
365  literature, which contains vector polynomials of degree up to (p+1).
366  For example, the RT_0 collection contains vector-valued linear functions
367  and, in particular, FiniteElementCollection::GetOrder() will,
368  correspondingly, return order 1. */
369  RT_FECollection(const int p, const int dim,
370  const int cb_type = BasisType::GaussLobatto,
371  const int ob_type = BasisType::GaussLegendre);
372 
374  Geometry::Type GeomType) const;
375  virtual int DofForGeometry(Geometry::Type GeomType) const
376  { return RT_dof[GeomType]; }
377  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
378  int Or) const;
379  virtual const char *Name() const { return rt_name; }
380  virtual int GetContType() const { return NORMAL; }
382 
383  int GetClosedBasisType() const { return cb_type; }
384  int GetOpenBasisType() const { return ob_type; }
385 
387  { return new RT_FECollection(p, dim, cb_type, ob_type); }
388 
389  virtual ~RT_FECollection();
390 };
391 
392 /** @brief Arbitrary order "H^{-1/2}-conforming" face finite elements defined on
393  the interface between mesh elements (faces); these are the normal trace FEs
394  of the H(div)-conforming FEs. */
396 {
397 public:
398  RT_Trace_FECollection(const int p, const int dim,
399  const int map_type = FiniteElement::INTEGRAL,
400  const int ob_type = BasisType::GaussLegendre);
401 };
402 
403 /** Arbitrary order discontinuous finite elements defined on the interface
404  between mesh elements (faces). The functions in this space are single-valued
405  on each face and are discontinuous across its boundary. */
407 {
408 public:
409  DG_Interface_FECollection(const int p, const int dim,
410  const int map_type = FiniteElement::VALUE,
411  const int ob_type = BasisType::GaussLegendre);
412 };
413 
414 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
416 {
417 protected:
418  int dim;
419  int cb_type; // closed BasisType
420  int ob_type; // open BasisType
421  char nd_name[32];
424  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
425 
426 public:
427  ND_FECollection(const int p, const int dim,
428  const int cb_type = BasisType::GaussLobatto,
429  const int ob_type = BasisType::GaussLegendre);
430 
431  virtual const FiniteElement *
433 
434  virtual int DofForGeometry(Geometry::Type GeomType) const
435  { return ND_dof[GeomType]; }
436 
437  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
438  int Or) const;
439  virtual const char *Name() const { return nd_name; }
440  virtual int GetContType() const { return TANGENTIAL; }
442 
443  int GetClosedBasisType() const { return cb_type; }
444  int GetOpenBasisType() const { return ob_type; }
445 
447  { return new ND_FECollection(p, dim, cb_type, ob_type); }
448 
449  virtual ~ND_FECollection();
450 };
451 
452 /** @brief Arbitrary order H(curl)-trace finite elements defined on the
453  interface between mesh elements (faces,edges); these are the tangential
454  trace FEs of the H(curl)-conforming FEs. */
456 {
457 public:
458  ND_Trace_FECollection(const int p, const int dim,
459  const int cb_type = BasisType::GaussLobatto,
460  const int ob_type = BasisType::GaussLegendre);
461 };
462 
463 /// Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
465 {
466 protected:
467  char nd_name[32];
470 
471 public:
472  ND_R1D_FECollection(const int p, const int dim,
473  const int cb_type = BasisType::GaussLobatto,
474  const int ob_type = BasisType::GaussLegendre);
475 
477  const
478  { return ND_Elements[GeomType]; }
479  virtual int DofForGeometry(Geometry::Type GeomType) const
480  { return ND_dof[GeomType]; }
481  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
482  int Or) const;
483  virtual const char *Name() const { return nd_name; }
484  virtual int GetContType() const { return TANGENTIAL; }
486 
487  virtual ~ND_R1D_FECollection();
488 };
489 
490 /// Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 1D.
492 {
493 protected:
494  char rt_name[32];
497 
498 public:
499  RT_R1D_FECollection(const int p, const int dim,
500  const int cb_type = BasisType::GaussLobatto,
501  const int ob_type = BasisType::GaussLegendre);
502 
504  const
505  { return RT_Elements[GeomType]; }
506  virtual int DofForGeometry(Geometry::Type GeomType) const
507  { return RT_dof[GeomType]; }
508  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
509  int Or) const;
510  virtual const char *Name() const { return rt_name; }
511  virtual int GetContType() const { return NORMAL; }
513 
514  virtual ~RT_R1D_FECollection();
515 };
516 
517 /// Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
519 {
520 protected:
521  char nd_name[32];
524  int *SegDofOrd[2];
525 
526 public:
527  ND_R2D_FECollection(const int p, const int dim,
528  const int cb_type = BasisType::GaussLobatto,
529  const int ob_type = BasisType::GaussLegendre);
530 
532  const
533  { return ND_Elements[GeomType]; }
534  virtual int DofForGeometry(Geometry::Type GeomType) const
535  { return ND_dof[GeomType]; }
536  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
537  int Or) const;
538  virtual const char *Name() const { return nd_name; }
539  virtual int GetContType() const { return TANGENTIAL; }
541 
542  virtual ~ND_R2D_FECollection();
543 };
544 
545 /** @brief Arbitrary order 3D H(curl)-trace finite elements in 2D defined on the
546  interface between mesh elements (edges); these are the tangential
547  trace FEs of the H(curl)-conforming FEs. */
549 {
550 public:
551  ND_R2D_Trace_FECollection(const int p, const int dim,
552  const int cb_type = BasisType::GaussLobatto,
553  const int ob_type = BasisType::GaussLegendre);
554 };
555 
556 /// Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 2D.
558 {
559 protected:
560  int ob_type; // open BasisType
561  char rt_name[32];
564  int *SegDofOrd[2];
565 
566  // Initialize only the face elements
567  void InitFaces(const int p, const int dim, const int map_type,
568  const bool signs);
569 
570  // Constructor used by the constructor of the RT_R2D_Trace_FECollection
571  RT_R2D_FECollection(const int p, const int dim, const int map_type,
572  const bool signs,
573  const int ob_type = BasisType::GaussLegendre);
574 
575 public:
576  RT_R2D_FECollection(const int p, const int dim,
577  const int cb_type = BasisType::GaussLobatto,
578  const int ob_type = BasisType::GaussLegendre);
579 
581  const
582  { return RT_Elements[GeomType]; }
583  virtual int DofForGeometry(Geometry::Type GeomType) const
584  { return RT_dof[GeomType]; }
585  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
586  int Or) const;
587  virtual const char *Name() const { return rt_name; }
588  virtual int GetContType() const { return NORMAL; }
590 
591  virtual ~RT_R2D_FECollection();
592 };
593 
594 /** @brief Arbitrary order 3D "H^{-1/2}-conforming" face finite elements defined on
595  the interface between mesh elements (faces); these are the normal trace FEs
596  of the H(div)-conforming FEs. */
598 {
599 public:
600  RT_R2D_Trace_FECollection(const int p, const int dim,
601  const int map_type = FiniteElement::INTEGRAL,
602  const int ob_type = BasisType::GaussLegendre);
603 };
604 
605 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
607 {
608 private:
609  NURBS1DFiniteElement *SegmentFE;
610  NURBS2DFiniteElement *QuadrilateralFE;
611  NURBS3DFiniteElement *ParallelepipedFE;
612 
613  mutable int mOrder; // >= 1 or VariableOrder
614  // The 'name' can be:
615  // 1) name = "NURBS" + "number", for fixed order, or
616  // 2) name = "NURBS", for VariableOrder.
617  // The name is updated before writing it to a stream, for example, see
618  // FiniteElementSpace::Save().
619  mutable char name[16];
620 
621 public:
622  enum { VariableOrder = -1 };
623 
624  /** @brief The parameter @a Order must be either a positive number, for fixed
625  order, or VariableOrder (default). */
626  explicit NURBSFECollection(int Order = VariableOrder);
627 
628  void Reset() const
629  {
630  SegmentFE->Reset();
631  QuadrilateralFE->Reset();
632  ParallelepipedFE->Reset();
633  }
634 
635  /** @brief Get the order of the NURBS collection: either a positive number,
636  when using fixed order, or VariableOrder. */
637  /** @note Not to be confused with FiniteElementCollection::GetOrder(). */
638  int GetOrder() const { return mOrder; }
639 
640  /** @brief Set the order and the name, based on the given @a Order: either a
641  positive number for fixed order, or VariableOrder. */
642  void SetOrder(int Order) const;
643 
644  virtual const FiniteElement *
646 
647  virtual int DofForGeometry(Geometry::Type GeomType) const;
648 
649  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
650  int Or) const;
651 
652  virtual const char *Name() const { return name; }
653 
654  virtual int GetContType() const { return CONTINUOUS; }
655 
657 
658  virtual ~NURBSFECollection();
659 };
660 
661 
662 /// Piecewise-(bi/tri)linear continuous finite elements.
664 {
665 private:
666  const PointFiniteElement PointFE;
667  const Linear1DFiniteElement SegmentFE;
668  const Linear2DFiniteElement TriangleFE;
669  const BiLinear2DFiniteElement QuadrilateralFE;
670  const Linear3DFiniteElement TetrahedronFE;
671  const TriLinear3DFiniteElement ParallelepipedFE;
672  const LinearWedgeFiniteElement WedgeFE;
673  const LinearPyramidFiniteElement PyramidFE;
674 public:
676 
677  virtual const FiniteElement *
679 
680  virtual int DofForGeometry(Geometry::Type GeomType) const;
681 
682  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
683  int Or) const;
684 
685  virtual const char * Name() const { return "Linear"; }
686 
687  virtual int GetContType() const { return CONTINUOUS; }
688 };
689 
690 /// Piecewise-(bi)quadratic continuous finite elements.
692 {
693 private:
694  const PointFiniteElement PointFE;
695  const Quad1DFiniteElement SegmentFE;
696  const Quad2DFiniteElement TriangleFE;
697  const BiQuad2DFiniteElement QuadrilateralFE;
698  const Quadratic3DFiniteElement TetrahedronFE;
699  const LagrangeHexFiniteElement ParallelepipedFE;
700  const H1_WedgeElement WedgeFE;
701 
702 public:
704  : FiniteElementCollection(2), ParallelepipedFE(2), WedgeFE(2) { }
705 
706  virtual const FiniteElement *
708 
709  virtual int DofForGeometry(Geometry::Type GeomType) const;
710 
711  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
712  int Or) const;
713 
714  virtual const char * Name() const { return "Quadratic"; }
715 
716  virtual int GetContType() const { return CONTINUOUS; }
717 };
718 
719 /// Version of QuadraticFECollection with positive basis functions.
721 {
722 private:
723  const QuadPos1DFiniteElement SegmentFE;
724  const BiQuadPos2DFiniteElement QuadrilateralFE;
725 
726 public:
728 
729  virtual const FiniteElement *
731 
732  virtual int DofForGeometry(Geometry::Type GeomType) const;
733 
734  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
735  int Or) const;
736 
737  virtual const char * Name() const { return "QuadraticPos"; }
738 
739  virtual int GetContType() const { return CONTINUOUS; }
740 };
741 
742 /// Piecewise-(bi)cubic continuous finite elements.
744 {
745 private:
746  const PointFiniteElement PointFE;
747  const Cubic1DFiniteElement SegmentFE;
748  const Cubic2DFiniteElement TriangleFE;
749  const BiCubic2DFiniteElement QuadrilateralFE;
750  const Cubic3DFiniteElement TetrahedronFE;
751  const LagrangeHexFiniteElement ParallelepipedFE;
752  const H1_WedgeElement WedgeFE;
753 
754 public:
757  ParallelepipedFE(3), WedgeFE(3, BasisType::ClosedUniform)
758  { }
759 
760  virtual const FiniteElement *
762 
763  virtual int DofForGeometry(Geometry::Type GeomType) const;
764 
765  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
766  int Or) const;
767 
768  virtual const char * Name() const { return "Cubic"; }
769 
770  virtual int GetContType() const { return CONTINUOUS; }
771 };
772 
773 /// Crouzeix-Raviart nonconforming elements in 2D.
775 {
776 private:
777  const P0SegmentFiniteElement SegmentFE;
778  const CrouzeixRaviartFiniteElement TriangleFE;
779  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
780 public:
782 
783  virtual const FiniteElement *
785 
786  virtual int DofForGeometry(Geometry::Type GeomType) const;
787 
788  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
789  int Or) const;
790 
791  virtual const char * Name() const { return "CrouzeixRaviart"; }
792 
793  virtual int GetContType() const { return DISCONTINUOUS; }
794 };
795 
796 /// Piecewise-linear nonconforming finite elements in 3D.
798 {
799 private:
800  const P0TriangleFiniteElement TriangleFE;
801  const P1TetNonConfFiniteElement TetrahedronFE;
802  const P0QuadFiniteElement QuadrilateralFE;
803  const RotTriLinearHexFiniteElement ParallelepipedFE;
804 
805 public:
807 
808  virtual const FiniteElement *
810 
811  virtual int DofForGeometry(Geometry::Type GeomType) const;
812 
813  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
814  int Or) const;
815 
816  virtual const char * Name() const { return "LinearNonConf3D"; }
817 
818  virtual int GetContType() const { return DISCONTINUOUS; }
819 };
820 
821 
822 /** @brief First order Raviart-Thomas finite elements in 2D. This class is kept
823  only for backward compatibility, consider using RT_FECollection instead. */
825 {
826 private:
827  const P0SegmentFiniteElement SegmentFE; // normal component on edge
828  const RT0TriangleFiniteElement TriangleFE;
829  const RT0QuadFiniteElement QuadrilateralFE;
830 public:
832 
833  virtual const FiniteElement *
835 
836  virtual int DofForGeometry(Geometry::Type GeomType) const;
837 
838  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
839  int Or) const;
840 
841  virtual const char * Name() const { return "RT0_2D"; }
842 
843  virtual int GetContType() const { return NORMAL; }
844 };
845 
846 /** @brief Second order Raviart-Thomas finite elements in 2D. This class is kept
847  only for backward compatibility, consider using RT_FECollection instead. */
849 {
850 private:
851  const P1SegmentFiniteElement SegmentFE; // normal component on edge
852  const RT1TriangleFiniteElement TriangleFE;
853  const RT1QuadFiniteElement QuadrilateralFE;
854 public:
856 
857  virtual const FiniteElement *
859 
860  virtual int DofForGeometry(Geometry::Type GeomType) const;
861 
862  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
863  int Or) const;
864 
865  virtual const char * Name() const { return "RT1_2D"; }
866 
867  virtual int GetContType() const { return NORMAL; }
868 };
869 
870 /** @brief Third order Raviart-Thomas finite elements in 2D. This class is kept
871  only for backward compatibility, consider using RT_FECollection instead. */
873 {
874 private:
875  const P2SegmentFiniteElement SegmentFE; // normal component on edge
876  const RT2TriangleFiniteElement TriangleFE;
877  const RT2QuadFiniteElement QuadrilateralFE;
878 public:
880 
881  virtual const FiniteElement *
883 
884  virtual int DofForGeometry(Geometry::Type GeomType) const;
885 
886  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
887  int Or) const;
888 
889  virtual const char * Name() const { return "RT2_2D"; }
890 
891  virtual int GetContType() const { return NORMAL; }
892 };
893 
894 /** @brief Piecewise-constant discontinuous finite elements in 2D. This class is
895  kept only for backward compatibility, consider using L2_FECollection
896  instead. */
898 {
899 private:
900  const P0TriangleFiniteElement TriangleFE;
901  const P0QuadFiniteElement QuadrilateralFE;
902 public:
904 
905  virtual const FiniteElement *
907 
908  virtual int DofForGeometry(Geometry::Type GeomType) const;
909 
910  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
911  int Or) const;
912 
913  virtual const char * Name() const { return "Const2D"; }
914 
915  virtual int GetContType() const { return DISCONTINUOUS; }
916 };
917 
918 /** @brief Piecewise-linear discontinuous finite elements in 2D. This class is
919  kept only for backward compatibility, consider using L2_FECollection
920  instead. */
922 {
923 private:
924  const Linear2DFiniteElement TriangleFE;
925  const BiLinear2DFiniteElement QuadrilateralFE;
926 
927 public:
929 
930  virtual const FiniteElement *
932 
933  virtual int DofForGeometry(Geometry::Type GeomType) const;
934 
935  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
936  int Or) const;
937 
938  virtual const char * Name() const { return "LinearDiscont2D"; }
939 
940  virtual int GetContType() const { return DISCONTINUOUS; }
941 };
942 
943 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
945 {
946 private:
947  // const CrouzeixRaviartFiniteElement TriangleFE;
948  const GaussLinear2DFiniteElement TriangleFE;
949  const GaussBiLinear2DFiniteElement QuadrilateralFE;
950 
951 public:
953 
954  virtual const FiniteElement *
956 
957  virtual int DofForGeometry(Geometry::Type GeomType) const;
958 
959  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
960  int Or) const;
961 
962  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
963 
964  virtual int GetContType() const { return DISCONTINUOUS; }
965 };
966 
967 /// Linear (P1) finite elements on quadrilaterals.
969 {
970 private:
971  const P1OnQuadFiniteElement QuadrilateralFE;
972 public:
974  virtual const FiniteElement *
976  virtual int DofForGeometry(Geometry::Type GeomType) const;
977  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
978  int Or) const;
979  virtual const char * Name() const { return "P1OnQuad"; }
980  virtual int GetContType() const { return DISCONTINUOUS; }
981 };
982 
983 /** @brief Piecewise-quadratic discontinuous finite elements in 2D. This class
984  is kept only for backward compatibility, consider using L2_FECollection
985  instead. */
987 {
988 private:
989  const Quad2DFiniteElement TriangleFE;
990  const BiQuad2DFiniteElement QuadrilateralFE;
991 
992 public:
994 
995  virtual const FiniteElement *
997 
998  virtual int DofForGeometry(Geometry::Type GeomType) const;
999 
1000  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1001  int Or) const;
1002 
1003  virtual const char * Name() const { return "QuadraticDiscont2D"; }
1004  virtual int GetContType() const { return DISCONTINUOUS; }
1005 };
1006 
1007 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
1009 {
1010 private:
1011  const BiQuadPos2DFiniteElement QuadrilateralFE;
1012 
1013 public:
1015  virtual const FiniteElement *
1016  FiniteElementForGeometry(Geometry::Type GeomType) const;
1017  virtual int DofForGeometry(Geometry::Type GeomType) const;
1018  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1019  int Or) const
1020  { return NULL; }
1021  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
1022  virtual int GetContType() const { return DISCONTINUOUS; }
1023 };
1024 
1025 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
1027 {
1028 private:
1029  // const Quad2DFiniteElement TriangleFE;
1030  const GaussQuad2DFiniteElement TriangleFE;
1031  const GaussBiQuad2DFiniteElement QuadrilateralFE;
1032 
1033 public:
1035 
1036  virtual const FiniteElement *
1037  FiniteElementForGeometry(Geometry::Type GeomType) const;
1038 
1039  virtual int DofForGeometry(Geometry::Type GeomType) const;
1040 
1041  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1042  int Or) const;
1043 
1044  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
1045  virtual int GetContType() const { return DISCONTINUOUS; }
1046 };
1047 
1048 /** @brief Piecewise-cubic discontinuous finite elements in 2D. This class is
1049  kept only for backward compatibility, consider using L2_FECollection
1050  instead. */
1052 {
1053 private:
1054  const Cubic2DFiniteElement TriangleFE;
1055  const BiCubic2DFiniteElement QuadrilateralFE;
1056 
1057 public:
1059 
1060  virtual const FiniteElement *
1061  FiniteElementForGeometry(Geometry::Type GeomType) const;
1062 
1063  virtual int DofForGeometry(Geometry::Type GeomType) const;
1064 
1065  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1066  int Or) const;
1067 
1068  virtual const char * Name() const { return "CubicDiscont2D"; }
1069  virtual int GetContType() const { return DISCONTINUOUS; }
1070 };
1071 
1072 /** @brief Piecewise-constant discontinuous finite elements in 3D. This class is
1073  kept only for backward compatibility, consider using L2_FECollection
1074  instead. */
1076 {
1077 private:
1078  const P0TetFiniteElement TetrahedronFE;
1079  const P0HexFiniteElement ParallelepipedFE;
1080  const P0WdgFiniteElement WedgeFE;
1081  const P0PyrFiniteElement PyramidFE;
1082 
1083 public:
1085 
1086  virtual const FiniteElement *
1087  FiniteElementForGeometry(Geometry::Type GeomType) const;
1088 
1089  virtual int DofForGeometry(Geometry::Type GeomType) const;
1090 
1091  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1092  int Or) const;
1093 
1094  virtual const char * Name() const { return "Const3D"; }
1095  virtual int GetContType() const { return DISCONTINUOUS; }
1096 };
1097 
1098 /** @brief Piecewise-linear discontinuous finite elements in 3D. This class is
1099  kept only for backward compatibility, consider using L2_FECollection
1100  instead. */
1102 {
1103 private:
1104  const Linear3DFiniteElement TetrahedronFE;
1105  const LinearPyramidFiniteElement PyramidFE;
1106  const LinearWedgeFiniteElement WedgeFE;
1107  const TriLinear3DFiniteElement ParallelepipedFE;
1108 
1109 public:
1111 
1112  virtual const FiniteElement *
1113  FiniteElementForGeometry(Geometry::Type GeomType) const;
1114 
1115  virtual int DofForGeometry(Geometry::Type GeomType) const;
1116 
1117  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1118  int Or) const;
1119 
1120  virtual const char * Name() const { return "LinearDiscont3D"; }
1121  virtual int GetContType() const { return DISCONTINUOUS; }
1122 };
1123 
1124 /** @brief Piecewise-quadratic discontinuous finite elements in 3D. This class
1125  is kept only for backward compatibility, consider using L2_FECollection
1126  instead. */
1128 {
1129 private:
1130  const Quadratic3DFiniteElement TetrahedronFE;
1131  const LagrangeHexFiniteElement ParallelepipedFE;
1132 
1133 public:
1135  : FiniteElementCollection(2), ParallelepipedFE(2) { }
1136 
1137  virtual const FiniteElement *
1138  FiniteElementForGeometry(Geometry::Type GeomType) const;
1139 
1140  virtual int DofForGeometry(Geometry::Type GeomType) const;
1141 
1142  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1143  int Or) const;
1144 
1145  virtual const char * Name() const { return "QuadraticDiscont3D"; }
1146  virtual int GetContType() const { return DISCONTINUOUS; }
1147 };
1148 
1149 /// Finite element collection on a macro-element.
1151 {
1152 private:
1153  const PointFiniteElement PointFE;
1154  const RefinedLinear1DFiniteElement SegmentFE;
1155  const RefinedLinear2DFiniteElement TriangleFE;
1156  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
1157  const RefinedLinear3DFiniteElement TetrahedronFE;
1158  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
1159 
1160 public:
1162 
1163  virtual const FiniteElement *
1164  FiniteElementForGeometry(Geometry::Type GeomType) const;
1165 
1166  virtual int DofForGeometry(Geometry::Type GeomType) const;
1167 
1168  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1169  int Or) const;
1170 
1171  virtual const char * Name() const { return "RefinedLinear"; }
1172  virtual int GetContType() const { return CONTINUOUS; }
1173 };
1174 
1175 /** @brief Lowest order Nedelec finite elements in 3D. This class is kept only
1176  for backward compatibility, consider using the new ND_FECollection
1177  instead. */
1179 {
1180 private:
1181  const Nedelec1HexFiniteElement HexahedronFE;
1182  const Nedelec1TetFiniteElement TetrahedronFE;
1183  const Nedelec1WdgFiniteElement WedgeFE;
1184  const Nedelec1PyrFiniteElement PyramidFE;
1185 
1186 public:
1188 
1189  virtual const FiniteElement *
1190  FiniteElementForGeometry(Geometry::Type GeomType) const;
1191 
1192  virtual int DofForGeometry(Geometry::Type GeomType) const;
1193 
1194  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1195  int Or) const;
1196 
1197  virtual const char * Name() const { return "ND1_3D"; }
1198  virtual int GetContType() const { return TANGENTIAL; }
1199 };
1200 
1201 /** @brief First order Raviart-Thomas finite elements in 3D. This class is kept
1202  only for backward compatibility, consider using RT_FECollection instead. */
1204 {
1205 private:
1206  const P0TriangleFiniteElement TriangleFE;
1207  const P0QuadFiniteElement QuadrilateralFE;
1208  const RT0HexFiniteElement HexahedronFE;
1209  const RT0TetFiniteElement TetrahedronFE;
1210  const RT0WdgFiniteElement WedgeFE;
1211  const RT0PyrFiniteElement PyramidFE;
1212 public:
1214 
1215  virtual const FiniteElement *
1216  FiniteElementForGeometry(Geometry::Type GeomType) const;
1217 
1218  virtual int DofForGeometry(Geometry::Type GeomType) const;
1219 
1220  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1221  int Or) const;
1222 
1223  virtual const char * Name() const { return "RT0_3D"; }
1224  virtual int GetContType() const { return NORMAL; }
1225 };
1226 
1227 /** @brief Second order Raviart-Thomas finite elements in 3D. This class is kept
1228  only for backward compatibility, consider using RT_FECollection instead. */
1230 {
1231 private:
1232  const Linear2DFiniteElement TriangleFE;
1233  const BiLinear2DFiniteElement QuadrilateralFE;
1234  const RT1HexFiniteElement HexahedronFE;
1235 public:
1237 
1238  virtual const FiniteElement *
1239  FiniteElementForGeometry(Geometry::Type GeomType) const;
1240 
1241  virtual int DofForGeometry(Geometry::Type GeomType) const;
1242 
1243  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1244  int Or) const;
1245 
1246  virtual const char * Name() const { return "RT1_3D"; }
1247  virtual int GetContType() const { return NORMAL; }
1248 };
1249 
1250 /// Discontinuous collection defined locally by a given finite element.
1252 {
1253 private:
1254  char d_name[32];
1255  Geometry::Type GeomType;
1256  FiniteElement *Local_Element;
1257 
1258 public:
1259  Local_FECollection(const char *fe_name);
1260 
1262  Geometry::Type GeomType_) const
1263  { return (GeomType == GeomType_) ? Local_Element : NULL; }
1264  virtual int DofForGeometry(Geometry::Type GeomType_) const
1265  { return (GeomType == GeomType_) ? Local_Element->GetDof() : 0; }
1266  virtual const int *DofOrderForOrientation(Geometry::Type GeomType_,
1267  int Or) const
1268  { return NULL; }
1269  virtual const char *Name() const { return d_name; }
1270 
1271  virtual ~Local_FECollection() { delete Local_Element; }
1272  virtual int GetContType() const { return DISCONTINUOUS; }
1273 };
1274 
1275 }
1276 
1277 #endif
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:869
Abstract class for all finite elements.
Definition: fe_base.hpp:235
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:638
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:986
A 2D 2nd order Raviart-Thomas vector element on a triangle.
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:606
virtual const char * Name() const
Definition: fe_coll.hpp:768
virtual const char * Name() const
Definition: fe_coll.hpp:1145
virtual const char * Name() const
Definition: fe_coll.hpp:1044
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1955
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1200
Piecewise-cubic discontinuous finite elements in 2D. This class is kept only for backward compatibili...
Definition: fe_coll.hpp:1051
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:1192
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:1026
int GetOpenBasisType() const
Definition: fe_coll.hpp:384
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:349
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1062
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1008
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:422
A 2D 1st order Raviart-Thomas vector element on a square.
const FiniteElement * GetFE(Geometry::Type geom, int p) const
Variable order version of FiniteElementForGeometry().
Definition: fe_coll.hpp:165
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:1018
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1241
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 FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2782
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:1157
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 const char * Name() const
Definition: fe_coll.hpp:1171
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:496
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:45
virtual const char * Name() const
Definition: fe_coll.hpp:1094
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2242
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:3376
int GetClosedBasisType() const
Definition: fe_coll.hpp:383
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:2256
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:2362
A 3D 1st order Nedelec element on a tetrahedron.
virtual const char * Name() const
Definition: fe_coll.hpp:587
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:605
virtual int DofForGeometry(Geometry::Type GeomType_) const
Definition: fe_coll.hpp:1264
virtual const char * Name() const
Definition: fe_coll.hpp:938
Field is discontinuous across element interfaces.
Definition: fe_coll.hpp:48
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:2511
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3410
static const int NumGeom
Definition: geom.hpp:42
A 3D 1st order Nedelec element on a pyramid.
virtual const char * Name() const
Definition: fe_coll.hpp:714
virtual int GetContType() const =0
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:2983
int GetClosedBasisType() const
Definition: fe_coll.hpp:443
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:464
virtual ~L2_FECollection()
Definition: fe_coll.cpp:2275
virtual const char * Name() const
Definition: fe_coll.hpp:979
virtual int GetContType() const
Definition: fe_coll.hpp:716
Normal component of vector field.
Definition: fe_coll.hpp:47
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1862
virtual const char * Name() const
Definition: fe_coll.hpp:510
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:3384
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:790
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:1266
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const =0
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:647
Piecewise-(bi/tri)linear continuous finite elements.
Definition: fe_coll.hpp:663
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:1876
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:959
void Reset() const
Definition: fe_nurbs.hpp:49
A 2D quadratic element on triangle with nodes at the vertices and midpoints of the triangle...
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:548
A 3D 0th order Raviert-Thomas element on a cube.
A 1D quadratic finite element with uniformly spaced nodes.
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1416
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
Definition: fe_coll.hpp:518
virtual const char * Name() const
Definition: fe_coll.hpp:379
virtual int GetContType() const
Definition: fe_coll.hpp:1172
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:476
virtual ~RT_R2D_FECollection()
Definition: fe_coll.cpp:3280
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2814
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:1024
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:693
const int base_p
Order as returned by GetOrder().
Definition: fe_coll.hpp:209
Second order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:848
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:236
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:2377
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2529
virtual const char * Name() const
Definition: fe_coll.hpp:652
A 3D constant element on a tetrahedron.
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:2578
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:3201
virtual const char * Name() const
Definition: fe_coll.hpp:841
A 3D Crouzeix-Raviart element on the tetrahedron.
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:446
A 2D bi-cubic element on a square with uniformly spaces nodes.
virtual const char * Name() const
Definition: fe_coll.hpp:1068
virtual const char * Name() const
Definition: fe_coll.hpp:538
Finite element collection on a macro-element.
Definition: fe_coll.hpp:1150
virtual int GetContType() const
Definition: fe_coll.hpp:1146
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:951
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
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition: fe_coll.hpp:199
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default)...
Definition: fe_coll.cpp:3352
virtual int GetContType() const
Definition: fe_coll.hpp:242
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1215
virtual int GetContType() const
Definition: fe_coll.hpp:1095
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:439
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1281
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:657
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:683
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:727
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2546
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1975
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:743
Class for finite elements with basis functions that return scalar values.
Definition: fe_base.hpp:629
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2989
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:855
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:580
virtual const char * Name() const
Definition: fe_coll.hpp:865
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1917
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:910
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:569
virtual int GetContType() const
Definition: fe_coll.hpp:1121
int HasFaceDofs(Geometry::Type geom, int p) const
Definition: fe_coll.cpp:25
virtual int GetContType() const
Definition: fe_coll.hpp:484
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:386
virtual int GetContType() const
Definition: fe_coll.hpp:380
virtual int GetContType() const
Definition: fe_coll.hpp:891
A 2D bi-linear element on a square with nodes at the &quot;Gaussian&quot; points.
A 1D refined linear element.
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:227
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:352
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:348
virtual int GetContType() const
Definition: fe_coll.hpp:940
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:937
A 2D constant element on a triangle.
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:322
int GetBasisType() const
Definition: fe_coll.hpp:329
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:468
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:583
virtual int GetContType() const
Definition: fe_coll.hpp:1198
A linear element defined on a triangular prism.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1402
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:1392
Piecewise-linear discontinuous finite elements in 2D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:921
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:522
int GetBasisType() const
Definition: fe_coll.hpp:243
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1551
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:2846
A 3D 0th order Raviert-Thomas element on a wedge.
Arbitrary order &quot;H^{1/2}-conforming&quot; trace finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:280
virtual int GetContType() const
Definition: fe_coll.hpp:440
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:3127
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2598
A 3D refined tri-linear element on a cube.
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:1251
Arbitrary order H(curl)-trace finite elements defined on the interface between mesh elements (faces...
Definition: fe_coll.hpp:455
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:3087
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:1008
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2834
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2497
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:2934
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:842
virtual int GetContType() const
Definition: fe_coll.hpp:867
Array< FiniteElementCollection * > var_orders
Definition: fe_coll.hpp:216
Third order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:872
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:994
First order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:1203
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:813
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:397
virtual int GetContType() const
Definition: fe_coll.hpp:588
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:1090
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1296
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:720
virtual ~Local_FECollection()
Definition: fe_coll.hpp:1271
virtual const char * Name() const
Definition: fe_coll.hpp:962
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:71
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:613
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1942
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:1352
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:1075
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:1314
A 3D linear element on a tetrahedron with nodes at the vertices of the tetrahedron.
A 1D linear element with nodes on the endpoints.
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:331
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:469
virtual int GetContType() const
Definition: fe_coll.hpp:321
A 2D Crouzeix-Raviart finite element on square.
A 2D linear element on triangle with nodes at the vertices of the triangle.
virtual int GetContType() const
Definition: fe_coll.hpp:1004
A quadratic element on triangle with nodes at the &quot;Gaussian&quot; points.
virtual const char * Name() const
Definition: fe_coll.hpp:241
virtual const char * Name() const
Definition: fe_coll.hpp:1120
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1898
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:375
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:897
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:226
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3264
const int * GetDofOrdering(Geometry::Type geom, int p, int ori) const
Variable order version of DofOrderForOrientation().
Definition: fe_coll.hpp:189
Piecewise-linear nonconforming finite elements in 3D.
Definition: fe_coll.hpp:797
Tangential components of vector field.
Definition: fe_coll.hpp:46
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:774
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:800
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:3097
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:562
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1335
int GetNumDof(Geometry::Type geom, int p) const
Variable order version of DofForGeometry().
Definition: fe_coll.hpp:177
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1141
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1032
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
Definition: fe_coll.hpp:338
H1Pos_FECollection(const int p, const int dim=3)
Definition: fe_coll.hpp:263
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:309
static void GetEdge(int &nv, v_t &v, int &ne, int &e, int &eo, const int edge_info)
Definition: fe_coll.cpp:332
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:323
H1Ser_FECollection(const int p, const int dim=2)
Definition: fe_coll.hpp:273
A 2D 2nd order Raviart-Thomas vector element on a square.
virtual const char * Name() const
Definition: fe_coll.hpp:791
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:1483
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const =0
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
virtual const char * Name() const
Definition: fe_coll.hpp:685
virtual const char * Name() const
Definition: fe_coll.hpp:889
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:3363
virtual const char * Name() const
Definition: fe_coll.hpp:816
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:896
virtual int GetContType() const
Definition: fe_coll.hpp:539
A 3D 1st order Raviert-Thomas element on a cube.
virtual const char * Name() const
Definition: fe_coll.hpp:1197
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
A 2D bi-quadratic element on a square with uniformly spaced 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:2796
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:2865
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
A 3D constant element on a wedge.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:503
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:776
Piecewise-quadratic discontinuous finite elements in 3D. This class is kept only for backward compati...
Definition: fe_coll.hpp:1127
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:630
virtual int GetContType() const
Definition: fe_coll.hpp:1224
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1322
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1449
A 3D constant element on a cube.
virtual const char * Name() const
Definition: fe_coll.hpp:1003
A 2D bi-linear element on a square with nodes at the vertices of the square.
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1376
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:944
virtual const char * Name() const
Definition: fe_coll.hpp:319
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:3004
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:51
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:534
A linear element on a triangle with nodes at the 3 &quot;Gaussian&quot; points.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1359
virtual int GetContType() const
Definition: fe_coll.hpp:793
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.cpp:296
virtual int GetContType() const
Definition: fe_coll.hpp:770
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:1178
Arbitrary order H1-conforming (continuous) finite elements with positive basis functions.
Definition: fe_coll.hpp:260
virtual const char * Name() const
Definition: fe_coll.hpp:483
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:1054
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:762
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1041
virtual const char * Name() const
Definition: fe_coll.hpp:65
virtual int GetContType() const
Definition: fe_coll.hpp:915
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1514
void InitVarOrder(int p) const
Definition: fe_coll.cpp:304
virtual int GetContType() const
Definition: fe_coll.hpp:739
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:434
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:252
virtual int GetContType() const
Definition: fe_coll.hpp:511
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:1529
A 2D bi-quadratic element on a square with nodes at the 9 &quot;Gaussian&quot; points.
virtual const char * Name() const
Definition: fe_coll.hpp:913
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:3403
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2920
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 int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:670
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:3215
Piecewise-linear discontinuous finite elements in 3D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:1101
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:3290
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1098
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1110
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:523
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1255
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:827
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:563
virtual const char * Name() const
Definition: fe_coll.hpp:1223
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:710
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1465
virtual int GetContType() const
Definition: fe_coll.hpp:654
virtual int GetContType() const
Definition: fe_coll.hpp:964
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:2558
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:3254
virtual int GetContType() const
Definition: fe_coll.hpp:1069
Tensor products of 1D Lagrange1DFiniteElement (only degree 2 is functional)
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:587
int GetOpenBasisType() const
Definition: fe_coll.hpp:444
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:691
First order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:824
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:1231
A 3D constant element on a pyramid.
A 2D constant element on a square.
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:506
virtual ~RT_R1D_FECollection()
Definition: fe_coll.cpp:2995
A linear element defined on a square pyramid.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:883
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:415
void Reset() const
Definition: fe_coll.hpp:628
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:479
A 2D cubic element on a triangle with uniformly spaced nodes.
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1075
A 1D cubic element with uniformly spaced nodes.
virtual const char * Name() const
Definition: fe_coll.hpp:1269
Arbitrary order &quot;H^{-1/2}-conforming&quot; face finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:395
virtual int GetContType() const
Definition: fe_coll.hpp:1272
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:557
virtual ~ND_R1D_FECollection()
Definition: fe_coll.cpp:2925
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
Second order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:1229
virtual const char * Name() const
Definition: fe_coll.hpp:1246
A 2D linear element on a square with 3 nodes at the vertices of the lower left triangle.
virtual const char * Name() const
Definition: fe_coll.hpp:1021
A 3D 0th order Raviert-Thomas element on a pyramid.
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType_) const
Definition: fe_coll.hpp:1261
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1126
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:3397
Field is continuous across element interfaces.
Definition: fe_coll.hpp:45
Piecewise-quadratic discontinuous finite elements in 2D. This class is kept only for backward compati...
Definition: fe_coll.hpp:986
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.hpp:531
virtual int GetContType() const
Definition: fe_coll.hpp:818
Arbitrary order 3D &quot;H^{-1/2}-conforming&quot; face finite elements defined on the interface between mesh e...
Definition: fe_coll.hpp:597
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:972
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:925
A 2D 3rd order Raviart-Thomas vector element on a triangle.
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:3310
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1178
virtual const char * Name() const
Definition: fe_coll.hpp:737
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:968
Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 1D.
Definition: fe_coll.hpp:491
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:1434
An arbitrary order 3D NURBS element on a cube.
Definition: fe_nurbs.hpp:115
virtual int GetContType() const
Definition: fe_coll.hpp:980
An arbitrary order 1D NURBS element on a segment.
Definition: fe_nurbs.hpp:63
virtual int GetContType() const
Definition: fe_coll.hpp:1247
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1165
virtual int GetContType() const
Definition: fe_coll.hpp:843
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:423
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:495
An arbitrary order 2D NURBS element on a square.
Definition: fe_nurbs.hpp:83
A 3D 1st order Nedelec element on a cube.
virtual int GetContType() const
Definition: fe_coll.hpp:687
A 2D Crouzeix-Raviart element on triangle.
virtual ~ND_R2D_FECollection()
Definition: fe_coll.cpp:3117
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1501
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:288
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:2914