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