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