MFEM  v4.3.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-2021, 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. */
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  { return H1_Elements[GeomType]; }
233  virtual int DofForGeometry(Geometry::Type GeomType) const
234  { return H1_dof[GeomType]; }
235  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
236  int Or) const;
237 
238  virtual const char *Name() const { return h1_name; }
239  virtual int GetContType() const { return CONTINUOUS; }
240  int GetBasisType() const { return b_type; }
241 
243 
244  /// Get the Cartesian to local H1 dof map
245  const int *GetDofMap(Geometry::Type GeomType) const;
246  /// Variable order version of GetDofMap
247  const int *GetDofMap(Geometry::Type GeomType, int p) const;
248 
250  { return new H1_FECollection(p, dim, b_type); }
251 
252  virtual ~H1_FECollection();
253 };
254 
255 /** @brief Arbitrary order H1-conforming (continuous) finite elements with
256  positive basis functions. */
258 {
259 public:
260  explicit H1Pos_FECollection(const int p, const int dim = 3)
261  : H1_FECollection(p, dim, BasisType::Positive) { }
262 };
263 
264 
265 /** Arbitrary order H1-conforming (continuous) serendipity finite elements;
266  Current implementation works in 2D only; 3D version is in development. */
268 {
269 public:
270  explicit H1Ser_FECollection(const int p, const int dim = 2)
271  : H1_FECollection(p, dim, BasisType::Serendipity) { };
272 };
273 
274 /** @brief Arbitrary order "H^{1/2}-conforming" trace finite elements defined on
275  the interface between mesh elements (faces,edges,vertices); these are the
276  trace FEs of the H1-conforming FEs. */
278 {
279 public:
280  H1_Trace_FECollection(const int p, const int dim,
281  const int btype = BasisType::GaussLobatto);
282 };
283 
284 /// Arbitrary order "L2-conforming" discontinuous finite elements.
286 {
287 private:
288  int dim;
289  int b_type; // BasisType
290  int m_type; // map type
291  char d_name[32];
294  int *SegDofOrd[2]; // for rotating segment dofs in 1D
295  int *TriDofOrd[6]; // for rotating triangle dofs in 2D
296  int *TetDofOrd[24]; // for rotating tetrahedron dofs in 3D
297  int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
298 
299 public:
300  L2_FECollection(const int p, const int dim,
301  const int btype = BasisType::GaussLegendre,
302  const int map_type = FiniteElement::VALUE);
303 
305  Geometry::Type GeomType) const
306  {
307  return L2_Elements[GeomType];
308  }
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  { return RT_Elements[GeomType]; }
376  virtual int DofForGeometry(Geometry::Type GeomType) const
377  { return RT_dof[GeomType]; }
378  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
379  int Or) const;
380  virtual const char *Name() const { return rt_name; }
381  virtual int GetContType() const { return NORMAL; }
383 
384  int GetClosedBasisType() const { return cb_type; }
385  int GetOpenBasisType() const { return ob_type; }
386 
388  { return new RT_FECollection(p, dim, cb_type, ob_type); }
389 
390  virtual ~RT_FECollection();
391 };
392 
393 /** @brief Arbitrary order "H^{-1/2}-conforming" face finite elements defined on
394  the interface between mesh elements (faces); these are the normal trace FEs
395  of the H(div)-conforming FEs. */
397 {
398 public:
399  RT_Trace_FECollection(const int p, const int dim,
400  const int map_type = FiniteElement::INTEGRAL,
401  const int ob_type = BasisType::GaussLegendre);
402 };
403 
404 /** Arbitrary order discontinuous finite elements defined on the interface
405  between mesh elements (faces). The functions in this space are single-valued
406  on each face and are discontinuous across its boundary. */
408 {
409 public:
410  DG_Interface_FECollection(const int p, const int dim,
411  const int map_type = FiniteElement::VALUE,
412  const int ob_type = BasisType::GaussLegendre);
413 };
414 
415 /// Arbitrary order H(curl)-conforming Nedelec finite elements.
417 {
418 protected:
419  int dim;
420  int cb_type; // closed BasisType
421  int ob_type; // open BasisType
422  char nd_name[32];
425  int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
426 
427 public:
428  ND_FECollection(const int p, const int dim,
429  const int cb_type = BasisType::GaussLobatto,
430  const int ob_type = BasisType::GaussLegendre);
431 
432  virtual const FiniteElement *
434  { return ND_Elements[GeomType]; }
435 
436  virtual int DofForGeometry(Geometry::Type GeomType) const
437  { return ND_dof[GeomType]; }
438 
439  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
440  int Or) const;
441  virtual const char *Name() const { return nd_name; }
442  virtual int GetContType() const { return TANGENTIAL; }
444 
445  int GetClosedBasisType() const { return cb_type; }
446  int GetOpenBasisType() const { return ob_type; }
447 
449  { return new ND_FECollection(p, dim, cb_type, ob_type); }
450 
451  virtual ~ND_FECollection();
452 };
453 
454 /** @brief Arbitrary order H(curl)-trace finite elements defined on the
455  interface between mesh elements (faces,edges); these are the tangential
456  trace FEs of the H(curl)-conforming FEs. */
458 {
459 public:
460  ND_Trace_FECollection(const int p, const int dim,
461  const int cb_type = BasisType::GaussLobatto,
462  const int ob_type = BasisType::GaussLegendre);
463 };
464 
465 /// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
467 {
468 private:
469  NURBS1DFiniteElement *SegmentFE;
470  NURBS2DFiniteElement *QuadrilateralFE;
471  NURBS3DFiniteElement *ParallelepipedFE;
472 
473  mutable int mOrder; // >= 1 or VariableOrder
474  // The 'name' can be:
475  // 1) name = "NURBS" + "number", for fixed order, or
476  // 2) name = "NURBS", for VariableOrder.
477  // The name is updated before writing it to a stream, for example, see
478  // FiniteElementSpace::Save().
479  mutable char name[16];
480 
481 public:
482  enum { VariableOrder = -1 };
483 
484  /** @brief The parameter @a Order must be either a positive number, for fixed
485  order, or VariableOrder (default). */
486  explicit NURBSFECollection(int Order = VariableOrder);
487 
488  void Reset() const
489  {
490  SegmentFE->Reset();
491  QuadrilateralFE->Reset();
492  ParallelepipedFE->Reset();
493  }
494 
495  /** @brief Get the order of the NURBS collection: either a positive number,
496  when using fixed order, or VariableOrder. */
497  /** @note Not to be confused with FiniteElementCollection::GetOrder(). */
498  int GetOrder() const { return mOrder; }
499 
500  /** @brief Set the order and the name, based on the given @a Order: either a
501  positive number for fixed order, or VariableOrder. */
502  void SetOrder(int Order) const;
503 
504  virtual const FiniteElement *
506 
507  virtual int DofForGeometry(Geometry::Type GeomType) const;
508 
509  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
510  int Or) const;
511 
512  virtual const char *Name() const { return name; }
513 
514  virtual int GetContType() const { return CONTINUOUS; }
515 
517 
518  virtual ~NURBSFECollection();
519 };
520 
521 
522 /// Piecewise-(bi/tri)linear continuous finite elements.
524 {
525 private:
526  const PointFiniteElement PointFE;
527  const Linear1DFiniteElement SegmentFE;
528  const Linear2DFiniteElement TriangleFE;
529  const BiLinear2DFiniteElement QuadrilateralFE;
530  const Linear3DFiniteElement TetrahedronFE;
531  const TriLinear3DFiniteElement ParallelepipedFE;
532  const H1_WedgeElement WedgeFE;
533 public:
535 
536  virtual const FiniteElement *
538 
539  virtual int DofForGeometry(Geometry::Type GeomType) const;
540 
541  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
542  int Or) const;
543 
544  virtual const char * Name() const { return "Linear"; }
545 
546  virtual int GetContType() const { return CONTINUOUS; }
547 };
548 
549 /// Piecewise-(bi)quadratic continuous finite elements.
551 {
552 private:
553  const PointFiniteElement PointFE;
554  const Quad1DFiniteElement SegmentFE;
555  const Quad2DFiniteElement TriangleFE;
556  const BiQuad2DFiniteElement QuadrilateralFE;
557  const Quadratic3DFiniteElement TetrahedronFE;
558  const LagrangeHexFiniteElement ParallelepipedFE;
559  const H1_WedgeElement WedgeFE;
560 
561 public:
563  : FiniteElementCollection(2), ParallelepipedFE(2), WedgeFE(2) { }
564 
565  virtual const FiniteElement *
567 
568  virtual int DofForGeometry(Geometry::Type GeomType) const;
569 
570  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
571  int Or) const;
572 
573  virtual const char * Name() const { return "Quadratic"; }
574 
575  virtual int GetContType() const { return CONTINUOUS; }
576 };
577 
578 /// Version of QuadraticFECollection with positive basis functions.
580 {
581 private:
582  const QuadPos1DFiniteElement SegmentFE;
583  const BiQuadPos2DFiniteElement QuadrilateralFE;
584 
585 public:
587 
588  virtual const FiniteElement *
590 
591  virtual int DofForGeometry(Geometry::Type GeomType) const;
592 
593  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
594  int Or) const;
595 
596  virtual const char * Name() const { return "QuadraticPos"; }
597 
598  virtual int GetContType() const { return CONTINUOUS; }
599 };
600 
601 /// Piecewise-(bi)cubic continuous finite elements.
603 {
604 private:
605  const PointFiniteElement PointFE;
606  const Cubic1DFiniteElement SegmentFE;
607  const Cubic2DFiniteElement TriangleFE;
608  const BiCubic2DFiniteElement QuadrilateralFE;
609  const Cubic3DFiniteElement TetrahedronFE;
610  const LagrangeHexFiniteElement ParallelepipedFE;
611  const H1_WedgeElement WedgeFE;
612 
613 public:
616  ParallelepipedFE(3), WedgeFE(3, BasisType::ClosedUniform)
617  { }
618 
619  virtual const FiniteElement *
621 
622  virtual int DofForGeometry(Geometry::Type GeomType) const;
623 
624  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
625  int Or) const;
626 
627  virtual const char * Name() const { return "Cubic"; }
628 
629  virtual int GetContType() const { return CONTINUOUS; }
630 };
631 
632 /// Crouzeix-Raviart nonconforming elements in 2D.
634 {
635 private:
636  const P0SegmentFiniteElement SegmentFE;
637  const CrouzeixRaviartFiniteElement TriangleFE;
638  const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
639 public:
641 
642  virtual const FiniteElement *
644 
645  virtual int DofForGeometry(Geometry::Type GeomType) const;
646 
647  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
648  int Or) const;
649 
650  virtual const char * Name() const { return "CrouzeixRaviart"; }
651 
652  virtual int GetContType() const { return DISCONTINUOUS; }
653 };
654 
655 /// Piecewise-linear nonconforming finite elements in 3D.
657 {
658 private:
659  const P0TriangleFiniteElement TriangleFE;
660  const P1TetNonConfFiniteElement TetrahedronFE;
661  const P0QuadFiniteElement QuadrilateralFE;
662  const RotTriLinearHexFiniteElement ParallelepipedFE;
663 
664 public:
666 
667  virtual const FiniteElement *
669 
670  virtual int DofForGeometry(Geometry::Type GeomType) const;
671 
672  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
673  int Or) const;
674 
675  virtual const char * Name() const { return "LinearNonConf3D"; }
676 
677  virtual int GetContType() const { return DISCONTINUOUS; }
678 };
679 
680 
681 /** @brief First order Raviart-Thomas finite elements in 2D. This class is kept
682  only for backward compatibility, consider using RT_FECollection instead. */
684 {
685 private:
686  const P0SegmentFiniteElement SegmentFE; // normal component on edge
687  const RT0TriangleFiniteElement TriangleFE;
688  const RT0QuadFiniteElement QuadrilateralFE;
689 public:
691 
692  virtual const FiniteElement *
694 
695  virtual int DofForGeometry(Geometry::Type GeomType) const;
696 
697  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
698  int Or) const;
699 
700  virtual const char * Name() const { return "RT0_2D"; }
701 
702  virtual int GetContType() const { return NORMAL; }
703 };
704 
705 /** @brief Second order Raviart-Thomas finite elements in 2D. This class is kept
706  only for backward compatibility, consider using RT_FECollection instead. */
708 {
709 private:
710  const P1SegmentFiniteElement SegmentFE; // normal component on edge
711  const RT1TriangleFiniteElement TriangleFE;
712  const RT1QuadFiniteElement QuadrilateralFE;
713 public:
715 
716  virtual const FiniteElement *
718 
719  virtual int DofForGeometry(Geometry::Type GeomType) const;
720 
721  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
722  int Or) const;
723 
724  virtual const char * Name() const { return "RT1_2D"; }
725 
726  virtual int GetContType() const { return NORMAL; }
727 };
728 
729 /** @brief Third order Raviart-Thomas finite elements in 2D. This class is kept
730  only for backward compatibility, consider using RT_FECollection instead. */
732 {
733 private:
734  const P2SegmentFiniteElement SegmentFE; // normal component on edge
735  const RT2TriangleFiniteElement TriangleFE;
736  const RT2QuadFiniteElement QuadrilateralFE;
737 public:
739 
740  virtual const FiniteElement *
742 
743  virtual int DofForGeometry(Geometry::Type GeomType) const;
744 
745  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
746  int Or) const;
747 
748  virtual const char * Name() const { return "RT2_2D"; }
749 
750  virtual int GetContType() const { return NORMAL; }
751 };
752 
753 /** @brief Piecewise-constant discontinuous finite elements in 2D. This class is
754  kept only for backward compatibility, consider using L2_FECollection
755  instead. */
757 {
758 private:
759  const P0TriangleFiniteElement TriangleFE;
760  const P0QuadFiniteElement QuadrilateralFE;
761 public:
763 
764  virtual const FiniteElement *
766 
767  virtual int DofForGeometry(Geometry::Type GeomType) const;
768 
769  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
770  int Or) const;
771 
772  virtual const char * Name() const { return "Const2D"; }
773 
774  virtual int GetContType() const { return DISCONTINUOUS; }
775 };
776 
777 /** @brief Piecewise-linear discontinuous finite elements in 2D. This class is
778  kept only for backward compatibility, consider using L2_FECollection
779  instead. */
781 {
782 private:
783  const Linear2DFiniteElement TriangleFE;
784  const BiLinear2DFiniteElement QuadrilateralFE;
785 
786 public:
788 
789  virtual const FiniteElement *
791 
792  virtual int DofForGeometry(Geometry::Type GeomType) const;
793 
794  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
795  int Or) const;
796 
797  virtual const char * Name() const { return "LinearDiscont2D"; }
798 
799  virtual int GetContType() const { return DISCONTINUOUS; }
800 };
801 
802 /// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
804 {
805 private:
806  // const CrouzeixRaviartFiniteElement TriangleFE;
807  const GaussLinear2DFiniteElement TriangleFE;
808  const GaussBiLinear2DFiniteElement QuadrilateralFE;
809 
810 public:
812 
813  virtual const FiniteElement *
815 
816  virtual int DofForGeometry(Geometry::Type GeomType) const;
817 
818  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
819  int Or) const;
820 
821  virtual const char * Name() const { return "GaussLinearDiscont2D"; }
822 
823  virtual int GetContType() const { return DISCONTINUOUS; }
824 };
825 
826 /// Linear (P1) finite elements on quadrilaterals.
828 {
829 private:
830  const P1OnQuadFiniteElement QuadrilateralFE;
831 public:
833  virtual const FiniteElement *
835  virtual int DofForGeometry(Geometry::Type GeomType) const;
836  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
837  int Or) const;
838  virtual const char * Name() const { return "P1OnQuad"; }
839  virtual int GetContType() const { return DISCONTINUOUS; }
840 };
841 
842 /** @brief Piecewise-quadratic discontinuous finite elements in 2D. This class
843  is kept only for backward compatibility, consider using L2_FECollection
844  instead. */
846 {
847 private:
848  const Quad2DFiniteElement TriangleFE;
849  const BiQuad2DFiniteElement QuadrilateralFE;
850 
851 public:
853 
854  virtual const FiniteElement *
856 
857  virtual int DofForGeometry(Geometry::Type GeomType) const;
858 
859  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
860  int Or) const;
861 
862  virtual const char * Name() const { return "QuadraticDiscont2D"; }
863  virtual int GetContType() const { return DISCONTINUOUS; }
864 };
865 
866 /// Version of QuadraticDiscont2DFECollection with positive basis functions.
868 {
869 private:
870  const BiQuadPos2DFiniteElement QuadrilateralFE;
871 
872 public:
874  virtual const FiniteElement *
876  virtual int DofForGeometry(Geometry::Type GeomType) const;
877  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
878  int Or) const
879  { return NULL; }
880  virtual const char * Name() const { return "QuadraticPosDiscont2D"; }
881  virtual int GetContType() const { return DISCONTINUOUS; }
882 };
883 
884 /// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
886 {
887 private:
888  // const Quad2DFiniteElement TriangleFE;
889  const GaussQuad2DFiniteElement TriangleFE;
890  const GaussBiQuad2DFiniteElement QuadrilateralFE;
891 
892 public:
894 
895  virtual const FiniteElement *
897 
898  virtual int DofForGeometry(Geometry::Type GeomType) const;
899 
900  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
901  int Or) const;
902 
903  virtual const char * Name() const { return "GaussQuadraticDiscont2D"; }
904  virtual int GetContType() const { return DISCONTINUOUS; }
905 };
906 
907 /** @brief Piecewise-cubic discontinuous finite elements in 2D. This class is
908  kept only for backward compatibility, consider using L2_FECollection
909  instead. */
911 {
912 private:
913  const Cubic2DFiniteElement TriangleFE;
914  const BiCubic2DFiniteElement QuadrilateralFE;
915 
916 public:
918 
919  virtual const FiniteElement *
921 
922  virtual int DofForGeometry(Geometry::Type GeomType) const;
923 
924  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
925  int Or) const;
926 
927  virtual const char * Name() const { return "CubicDiscont2D"; }
928  virtual int GetContType() const { return DISCONTINUOUS; }
929 };
930 
931 /** @brief Piecewise-constant discontinuous finite elements in 3D. This class is
932  kept only for backward compatibility, consider using L2_FECollection
933  instead. */
935 {
936 private:
937  const P0TetFiniteElement TetrahedronFE;
938  const P0HexFiniteElement ParallelepipedFE;
939  const L2_WedgeElement WedgeFE;
940 
941 public:
943 
944  virtual const FiniteElement *
946 
947  virtual int DofForGeometry(Geometry::Type GeomType) const;
948 
949  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
950  int Or) const;
951 
952  virtual const char * Name() const { return "Const3D"; }
953  virtual int GetContType() const { return DISCONTINUOUS; }
954 };
955 
956 /** @brief Piecewise-linear discontinuous finite elements in 3D. This class is
957  kept only for backward compatibility, consider using L2_FECollection
958  instead. */
960 {
961 private:
962  const Linear3DFiniteElement TetrahedronFE;
963  const TriLinear3DFiniteElement ParallelepipedFE;
964 
965 public:
967 
968  virtual const FiniteElement *
970 
971  virtual int DofForGeometry(Geometry::Type GeomType) const;
972 
973  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
974  int Or) const;
975 
976  virtual const char * Name() const { return "LinearDiscont3D"; }
977  virtual int GetContType() const { return DISCONTINUOUS; }
978 };
979 
980 /** @brief Piecewise-quadratic discontinuous finite elements in 3D. This class
981  is kept only for backward compatibility, consider using L2_FECollection
982  instead. */
984 {
985 private:
986  const Quadratic3DFiniteElement TetrahedronFE;
987  const LagrangeHexFiniteElement ParallelepipedFE;
988 
989 public:
991  : FiniteElementCollection(2), ParallelepipedFE(2) { }
992 
993  virtual const FiniteElement *
995 
996  virtual int DofForGeometry(Geometry::Type GeomType) const;
997 
998  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
999  int Or) const;
1000 
1001  virtual const char * Name() const { return "QuadraticDiscont3D"; }
1002  virtual int GetContType() const { return DISCONTINUOUS; }
1003 };
1004 
1005 /// Finite element collection on a macro-element.
1007 {
1008 private:
1009  const PointFiniteElement PointFE;
1010  const RefinedLinear1DFiniteElement SegmentFE;
1011  const RefinedLinear2DFiniteElement TriangleFE;
1012  const RefinedBiLinear2DFiniteElement QuadrilateralFE;
1013  const RefinedLinear3DFiniteElement TetrahedronFE;
1014  const RefinedTriLinear3DFiniteElement ParallelepipedFE;
1015 
1016 public:
1018 
1019  virtual const FiniteElement *
1020  FiniteElementForGeometry(Geometry::Type GeomType) const;
1021 
1022  virtual int DofForGeometry(Geometry::Type GeomType) const;
1023 
1024  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1025  int Or) const;
1026 
1027  virtual const char * Name() const { return "RefinedLinear"; }
1028  virtual int GetContType() const { return CONTINUOUS; }
1029 };
1030 
1031 /** @brief Lowest order Nedelec finite elements in 3D. This class is kept only
1032  for backward compatibility, consider using the new ND_FECollection
1033  instead. */
1035 {
1036 private:
1037  const Nedelec1HexFiniteElement HexahedronFE;
1038  const Nedelec1TetFiniteElement TetrahedronFE;
1039 
1040 public:
1042 
1043  virtual const FiniteElement *
1044  FiniteElementForGeometry(Geometry::Type GeomType) const;
1045 
1046  virtual int DofForGeometry(Geometry::Type GeomType) const;
1047 
1048  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1049  int Or) const;
1050 
1051  virtual const char * Name() const { return "ND1_3D"; }
1052  virtual int GetContType() const { return TANGENTIAL; }
1053 };
1054 
1055 /** @brief First order Raviart-Thomas finite elements in 3D. This class is kept
1056  only for backward compatibility, consider using RT_FECollection instead. */
1058 {
1059 private:
1060  const P0TriangleFiniteElement TriangleFE;
1061  const P0QuadFiniteElement QuadrilateralFE;
1062  const RT0HexFiniteElement HexahedronFE;
1063  const RT0TetFiniteElement TetrahedronFE;
1064 public:
1066 
1067  virtual const FiniteElement *
1068  FiniteElementForGeometry(Geometry::Type GeomType) const;
1069 
1070  virtual int DofForGeometry(Geometry::Type GeomType) const;
1071 
1072  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1073  int Or) const;
1074 
1075  virtual const char * Name() const { return "RT0_3D"; }
1076  virtual int GetContType() const { return NORMAL; }
1077 };
1078 
1079 /** @brief Second order Raviart-Thomas finite elements in 3D. This class is kept
1080  only for backward compatibility, consider using RT_FECollection instead. */
1082 {
1083 private:
1084  const Linear2DFiniteElement TriangleFE;
1085  const BiLinear2DFiniteElement QuadrilateralFE;
1086  const RT1HexFiniteElement HexahedronFE;
1087 public:
1089 
1090  virtual const FiniteElement *
1091  FiniteElementForGeometry(Geometry::Type GeomType) const;
1092 
1093  virtual int DofForGeometry(Geometry::Type GeomType) const;
1094 
1095  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1096  int Or) const;
1097 
1098  virtual const char * Name() const { return "RT1_3D"; }
1099  virtual int GetContType() const { return NORMAL; }
1100 };
1101 
1102 /// Discontinuous collection defined locally by a given finite element.
1104 {
1105 private:
1106  char d_name[32];
1107  Geometry::Type GeomType;
1108  FiniteElement *Local_Element;
1109 
1110 public:
1111  Local_FECollection(const char *fe_name);
1112 
1114  Geometry::Type GeomType_) const
1115  { return (GeomType == GeomType_) ? Local_Element : NULL; }
1116  virtual int DofForGeometry(Geometry::Type GeomType_) const
1117  { return (GeomType == GeomType_) ? Local_Element->GetDof() : 0; }
1118  virtual const int *DofOrderForOrientation(Geometry::Type GeomType,
1119  int Or) const
1120  { return NULL; }
1121  virtual const char *Name() const { return d_name; }
1122 
1123  virtual ~Local_FECollection() { delete Local_Element; }
1124  virtual int GetContType() const { return DISCONTINUOUS; }
1125 };
1126 
1127 }
1128 
1129 #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:864
Abstract class for all finite elements.
Definition: fe.hpp:243
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order...
Definition: fe_coll.hpp:498
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:981
A 2D 2nd order Raviart-Thomas vector element on a triangle.
Definition: fe.hpp:1510
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:466
Field is discontinuous across element interfaces.
Definition: fe_coll.hpp:48
virtual const char * Name() const
Definition: fe_coll.hpp:627
virtual const char * Name() const
Definition: fe_coll.hpp:1001
virtual const char * Name() const
Definition: fe_coll.hpp:903
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1919
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1195
Arbitrary order L2 elements in 3D on a wedge.
Definition: fe.hpp:2763
Piecewise-cubic discontinuous finite elements in 2D. This class is kept only for backward compatibili...
Definition: fe_coll.hpp:910
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:1187
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:885
int GetOpenBasisType() const
Definition: fe_coll.hpp:385
int RT_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:349
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1057
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1003
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:423
A 2D 1st order Raviart-Thomas vector element on a square.
Definition: fe.hpp:1481
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:877
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1236
A 1D constant element on a segment.
Definition: fe.hpp:1441
A 1D linear element with nodes at 1/3 and 2/3 (trace of RT1)
Definition: fe.hpp:1619
A 1D quadratic element with nodes at the Gaussian points (trace of RT2)
Definition: fe.hpp:1630
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:1152
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:1266
virtual const char * Name() const
Definition: fe_coll.hpp:1027
virtual FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:42
virtual const char * Name() const
Definition: fe_coll.hpp:952
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:230
virtual ~NURBSFECollection()
Definition: fe_coll.cpp:2826
int GetClosedBasisType() const
Definition: fe_coll.hpp:384
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:2202
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:2302
A 3D 1st order Nedelec element on a tetrahedron.
Definition: fe.hpp:1835
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:600
virtual int DofForGeometry(Geometry::Type GeomType_) const
Definition: fe_coll.hpp:1116
virtual const char * Name() const
Definition: fe_coll.hpp:797
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:2436
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2860
static const int NumGeom
Definition: geom.hpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:573
virtual int GetContType() const =0
Tangential components of vector field.
Definition: fe_coll.hpp:46
int GetClosedBasisType() const
Definition: fe_coll.hpp:445
const Geometry::Type geom
Definition: ex1.cpp:40
A 2D 3rd order Raviart-Thomas vector element on a square.
Definition: fe.hpp:1588
virtual ~L2_FECollection()
Definition: fe_coll.cpp:2221
virtual const char * Name() const
Definition: fe_coll.hpp:838
virtual int GetContType() const
Definition: fe_coll.hpp:575
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2834
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:785
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:642
Piecewise-(bi/tri)linear continuous finite elements.
Definition: fe_coll.hpp:523
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:1840
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:954
void Reset() const
Definition: fe.hpp:3392
A 2D quadratic element on triangle with nodes at the vertices and midpoints of the triangle...
Definition: fe.hpp:1162
A 0D point finite element.
Definition: fe.hpp:1010
A 3D 0th order Raviert-Thomas element on a cube.
Definition: fe.hpp:1859
A 1D quadratic finite element with uniformly spaced nodes.
Definition: fe.hpp:1130
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1403
virtual const char * Name() const
Definition: fe_coll.hpp:380
virtual int GetContType() const
Definition: fe_coll.hpp:1028
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2709
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:1019
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:688
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:707
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:233
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition: fe_coll.cpp:2317
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:2454
virtual const char * Name() const
Definition: fe_coll.hpp:512
A 3D constant element on a tetrahedron.
Definition: fe.hpp:1668
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:2502
virtual const char * Name() const
Definition: fe_coll.hpp:700
A 3D Crouzeix-Raviart element on the tetrahedron.
Definition: fe.hpp:1657
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:448
A 2D bi-cubic element on a square with uniformly spaces nodes.
Definition: fe.hpp:1259
virtual const char * Name() const
Definition: fe_coll.hpp:927
Finite element collection on a macro-element.
Definition: fe_coll.hpp:1006
virtual int GetContType() const
Definition: fe_coll.hpp:1002
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:946
A 2D refined bi-linear FE on a square.
Definition: fe.hpp:1770
Possible basis types. Note that not all elements can use all BasisType(s).
Definition: fe.hpp:27
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:2802
virtual int GetContType() const
Definition: fe_coll.hpp:239
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1210
virtual int GetContType() const
Definition: fe_coll.hpp:953
virtual int DofForGeometry(Geometry::Type GeomType) const =0
virtual const char * Name() const
Definition: fe_coll.hpp:441
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1274
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:652
virtual const int * DofOrderForOrientation(Geometry::Type GeomType, int Or) const
Returns an array, say p, that maps a local permuted index i to a local base index: base_i = p[i]...
Definition: fe_coll.cpp:678
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:722
virtual ~RT_FECollection()
Definition: fe_coll.cpp:2470
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE)
Definition: fe_coll.cpp:1939
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:304
Piecewise-(bi)cubic continuous finite elements.
Definition: fe_coll.hpp:602
Class for finite elements with basis functions that return scalar values.
Definition: fe.hpp:629
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:850
virtual const char * Name() const
Definition: fe_coll.hpp:724
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition: fe_coll.cpp:1881
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:905
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:566
virtual int GetContType() const
Definition: fe_coll.hpp:977
int HasFaceDofs(Geometry::Type geom, int p) const
Definition: fe_coll.cpp:25
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:387
virtual int GetContType() const
Definition: fe_coll.hpp:381
virtual int GetContType() const
Definition: fe_coll.hpp:750
A 2D bi-linear element on a square with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:1101
A 1D refined linear element.
Definition: fe.hpp:1717
int H1_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:223
A 2D refined linear element on a triangle.
Definition: fe.hpp:1737
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
A 2D refined linear element on a tetrahedron.
Definition: fe.hpp:1757
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:349
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:348
virtual int GetContType() const
Definition: fe_coll.hpp:799
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:932
A 2D constant element on a triangle.
Definition: fe.hpp:1317
static void GetNVE(int &nv, int &ne)
Definition: fe_coll.cpp:319
int GetBasisType() const
Definition: fe_coll.hpp:329
virtual int GetContType() const
Definition: fe_coll.hpp:1052
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1391
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:1381
Piecewise-linear discontinuous finite elements in 2D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:780
int GetBasisType() const
Definition: fe_coll.hpp:240
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto)
Definition: fe_coll.cpp:1532
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:2741
Arbitrary order &quot;H^{1/2}-conforming&quot; trace finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:277
virtual int GetContType() const
Definition: fe_coll.hpp:442
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition: fe_coll.cpp:2522
A 3D refined tri-linear element on a cube.
Definition: fe.hpp:1790
Discontinuous collection defined locally by a given finite element.
Definition: fe_coll.hpp:1103
Arbitrary order H(curl)-trace finite elements defined on the interface between mesh elements (faces...
Definition: fe_coll.hpp:457
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition: fe_coll.hpp:867
virtual ~ND_FECollection()
Definition: fe_coll.cpp:2729
A 3D quadratic element on a tetrahedron with uniformly spaced nodes.
Definition: fe.hpp:1379
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:837
virtual int GetContType() const
Definition: fe_coll.hpp:726
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:731
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:989
First order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp:1057
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:808
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:394
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:1085
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1287
Version of QuadraticFECollection with positive basis functions.
Definition: fe_coll.hpp:579
virtual ~Local_FECollection()
Definition: fe_coll.hpp:1123
virtual const char * Name() const
Definition: fe_coll.hpp:821
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:608
virtual ~H1_FECollection()
Definition: fe_coll.cpp:1906
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:1341
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
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:1303
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:1118
A 3D linear element on a tetrahedron with nodes at the vertices of the tetrahedron.
Definition: fe.hpp:1350
A 1D linear element with nodes on the endpoints.
Definition: fe.hpp:1023
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:331
virtual int GetContType() const
Definition: fe_coll.hpp:321
A 2D Crouzeix-Raviart finite element on square.
Definition: fe.hpp:1429
A 2D linear element on triangle with nodes at the vertices of the triangle.
Definition: fe.hpp:1043
virtual int GetContType() const
Definition: fe_coll.hpp:863
A quadratic element on triangle with nodes at the &quot;Gaussian&quot; points.
Definition: fe.hpp:1186
virtual const char * Name() const
Definition: fe_coll.hpp:238
virtual const char * Name() const
Definition: fe_coll.hpp:976
FiniteElementCollection * GetTraceCollection() const
Definition: fe_coll.cpp:1862
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:376
Piecewise-constant discontinuous finite elements in 2D. This class is kept only for backward compatib...
Definition: fe_coll.hpp:756
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition: fe_coll.hpp:222
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:656
Crouzeix-Raviart nonconforming elements in 2D.
Definition: fe_coll.hpp:633
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:795
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1324
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:1136
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:1027
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:260
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:329
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:270
A 2D 2nd order Raviart-Thomas vector element on a square.
Definition: fe.hpp:1539
virtual const char * Name() const
Definition: fe_coll.hpp:650
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:1464
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:544
virtual const char * Name() const
Definition: fe_coll.hpp:748
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:2813
virtual const char * Name() const
Definition: fe_coll.hpp:675
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:891
A 3D 1st order Raviert-Thomas element on a cube.
Definition: fe.hpp:1889
virtual const char * Name() const
Definition: fe_coll.hpp:1051
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:323
A 2D bi-quadratic element on a square with uniformly spaced nodes.
Definition: fe.hpp:1203
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:2686
A 3D tri-linear element on a cube with nodes at the vertices of the cube.
Definition: fe.hpp:1392
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:771
Piecewise-quadratic discontinuous finite elements in 3D. This class is kept only for backward compati...
Definition: fe_coll.hpp:983
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:625
Field is continuous across element interfaces.
Definition: fe_coll.hpp:45
virtual int GetContType() const
Definition: fe_coll.hpp:1076
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1311
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1434
A 3D constant element on a cube.
Definition: fe.hpp:1681
virtual const char * Name() const
Definition: fe_coll.hpp:862
A 2D bi-linear element on a square with nodes at the vertices of the square.
Definition: fe.hpp:1065
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1365
Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
Definition: fe_coll.hpp:803
virtual const char * Name() const
Definition: fe_coll.hpp:319
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
Definition: fe_coll.cpp:48
A linear element on a triangle with nodes at the 3 &quot;Gaussian&quot; points.
Definition: fe.hpp:1089
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1348
virtual int GetContType() const
Definition: fe_coll.hpp:652
virtual FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.cpp:293
virtual int GetContType() const
Definition: fe_coll.hpp:629
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:1034
Arbitrary order H1-conforming (continuous) finite elements with positive basis functions.
Definition: fe_coll.hpp:257
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:1049
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:757
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1036
virtual const char * Name() const
Definition: fe_coll.hpp:61
virtual int GetContType() const
Definition: fe_coll.hpp:774
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1495
void InitVarOrder(int p) const
Definition: fe_coll.cpp:301
virtual int GetContType() const
Definition: fe_coll.hpp:598
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:436
FiniteElementCollection * Clone(int p) const
Instantiate a new collection of the same type with a different order.
Definition: fe_coll.hpp:249
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:1510
A 2D bi-quadratic element on a square with nodes at the 9 &quot;Gaussian&quot; points.
Definition: fe.hpp:1246
virtual const char * Name() const
Definition: fe_coll.hpp:772
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:2853
A 1D quadratic positive element utilizing the 2nd order Bernstein basis.
Definition: fe.hpp:1150
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:665
Piecewise-linear discontinuous finite elements in 3D. This class is kept only for backward compatibil...
Definition: fe_coll.hpp:959
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1093
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1105
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:373
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1249
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:822
virtual const char * Name() const
Definition: fe_coll.hpp:1075
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:705
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1448
virtual int GetContType() const
Definition: fe_coll.hpp:514
virtual int GetContType() const
Definition: fe_coll.hpp:823
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:2482
virtual int GetContType() const
Definition: fe_coll.hpp:928
Tensor products of 1D Lagrange1DFiniteElement (only degree 2 is functional)
Definition: fe.hpp:1695
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:583
int GetOpenBasisType() const
Definition: fe_coll.hpp:446
Piecewise-(bi)quadratic continuous finite elements.
Definition: fe_coll.hpp:550
First order Raviart-Thomas finite elements in 2D. This class is kept only for backward compatibility...
Definition: fe_coll.hpp: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:1226
A 2D constant element on a square.
Definition: fe.hpp:1335
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:878
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:416
void Reset() const
Definition: fe_coll.hpp:488
A 2D cubic element on a triangle with uniformly spaced nodes.
Definition: fe.hpp:1287
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1070
A 1D cubic element with uniformly spaced nodes.
Definition: fe.hpp:1274
virtual const char * Name() const
Definition: fe_coll.hpp:1121
Arbitrary order &quot;H^{-1/2}-conforming&quot; face finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:396
virtual int GetContType() const
Definition: fe_coll.hpp:1124
A 2D 1st order Raviart-Thomas vector element on a triangle.
Definition: fe.hpp:1452
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:1081
virtual const char * Name() const
Definition: fe_coll.hpp:1098
A 2D linear element on a square with 3 nodes at the vertices of the lower left triangle.
Definition: fe.hpp:1117
virtual const char * Name() const
Definition: fe_coll.hpp:880
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType_) const
Definition: fe_coll.hpp:1113
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1121
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:2847
Piecewise-quadratic discontinuous finite elements in 2D. This class is kept only for backward compati...
Definition: fe_coll.hpp:845
Arbitrary order H1 elements in 3D on a wedge.
Definition: fe.hpp:2500
A 3D 0th order Raviert-Thomas element on a tetrahedron.
Definition: fe.hpp:1919
virtual int GetContType() const
Definition: fe_coll.hpp:677
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:967
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:920
A 2D 3rd order Raviart-Thomas vector element on a triangle.
Definition: fe.hpp:1568
Local_FECollection(const char *fe_name)
Definition: fe_coll.cpp:2760
virtual int DofForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1173
virtual const char * Name() const
Definition: fe_coll.hpp:596
Linear (P1) finite elements on quadrilaterals.
Definition: fe_coll.hpp:827
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:1419
An arbitrary order 3D NURBS element on a cube.
Definition: fe.hpp:3458
virtual int GetContType() const
Definition: fe_coll.hpp:839
An arbitrary order 1D NURBS element on a segment.
Definition: fe.hpp:3406
virtual int GetContType() const
Definition: fe_coll.hpp:1099
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1160
virtual int GetContType() const
Definition: fe_coll.hpp:702
int ND_dof[Geometry::NumGeom]
Definition: fe_coll.hpp:424
An arbitrary order 2D NURBS element on a square.
Definition: fe.hpp:3426
A 3D 1st order Nedelec element on a cube.
Definition: fe.hpp:1811
virtual int GetContType() const
Definition: fe_coll.hpp:546
A 2D Crouzeix-Raviart element on triangle.
Definition: fe.hpp:1416
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.cpp:1482
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const
Definition: fe_coll.hpp:433