MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
fe_coll.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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 | ND_R1D_[DIM]_[ORDER] | H(curl) | * | 1 / 0 | H_CURL | 3D H(curl)-conforming Nedelec vector elements in 1D. |
123 | ND_R1D@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(curl) | * | * / * | H_CURL | 3D H(curl)-conforming Nedelec vector elements in 1D. |
124 | ND_R2D_[DIM]_[ORDER] | H(curl) | * | 1 / 0 | H_CURL | 3D H(curl)-conforming Nedelec vector elements in 2D. |
125 | ND_R2D@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(curl) | * | * / * | H_CURL | 3D H(curl)-conforming Nedelec vector elements in 2D. |
126 | RT_[DIM]_[ORDER] | H(div) | * | 1 / 0 | H_DIV | Raviart-Thomas vector elements |
127 | RT@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(div) | * | * / * | H_DIV | Raviart-Thomas vector elements |
128 | 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) |
129 | 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) |
130 | 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) |
131 | 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) |
132 | RT_R1D_[DIM]_[ORDER] | H(div) | * | 1 / 0 | H_DIV | 3D H(div)-conforming Raviart-Thomas vector elements in 1D. |
133 | RT_R1D@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(div) | * | * / * | H_DIV | 3D H(div)-conforming Raviart-Thomas vector elements in 1D. |
134 | RT_R2D_[DIM]_[ORDER] | H(div) | * | 1 / 0 | H_DIV | 3D H(div)-conforming Raviart-Thomas vector elements in 2D. |
135 | RT_R2D@[CBTYPE][OBTYPE]_[DIM]_[ORDER] | H(div) | * | * / * | H_DIV | 3D H(div)-conforming Raviart-Thomas vector elements in 2D. |
136 | L2_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinuous L2 elements |
137 | L2_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | VALUE | Discontinuous L2 elements |
138 | L2Int_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinuous L2 elements |
139 | L2Int_T[BTYPE]_[DIM]_[ORDER] | L2 | * | 0 | INTEGRAL | Discontinuous L2 elements |
140 | DG_Iface_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
141 | DG_Iface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | VALUE | Discontinuous elements on the interface between mesh elements (faces) |
142 | DG_IntIface_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
143 | DG_IntIface@[BTYPE]_[DIM]_[ORDER] | - | * | 0 | INTEGRAL | Discontinuous elements on the interface between mesh elements (faces) |
144 | NURBS[ORDER] | - | * | - | VALUE | Non-Uniform Rational B-Splines (NURBS) elements |
145 | LinearNonConf3D | - | 1 | 1 | VALUE | Piecewise-linear nonconforming finite elements in 3D |
146 | CrouzeixRaviart | - | - | - | - | Crouzeix-Raviart nonconforming elements in 2D |
147 | Local_[FENAME] | - | - | - | - | Special collection that builds a local version out of the FENAME collection |
148 |-|-|-|-|-|-|
149 | Linear | H1 | 1 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
150 | Quadratic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
151 | QuadraticPos | H1 | 2 | 2 | VALUE | Left in for backward compatibility, consider using H1_ |
152 | Cubic | H1 | 2 | 1 | VALUE | Left in for backward compatibility, consider using H1_ |
153 | Const2D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
154 | Const3D | L2 | 0 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
155 | LinearDiscont2D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
156 | GaussLinearDiscont2D | L2 | 1 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
157 | P1OnQuad | H1 | 1 | 1 | VALUE | Linear P1 element with 3 nodes on a square |
158 | QuadraticDiscont2D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
159 | QuadraticPosDiscont2D | L2 | 2 | 2 | VALUE | Left in for backward compatibility, consider using L2_ |
160 | GaussQuadraticDiscont2D | L2 | 2 | 0 | VALUE | Left in for backward compatibility, consider using L2_ |
161 | CubicDiscont2D | L2 | 3 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
162 | LinearDiscont3D | L2 | 1 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
163 | QuadraticDiscont3D | L2 | 2 | 1 | VALUE | Left in for backward compatibility, consider using L2_ |
164 | ND1_3D | H(Curl) | 1 | 1 / 0 | H_CURL | Left in for backward compatibility, consider using ND_ |
165 | RT0_2D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
166 | RT1_2D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
167 | RT2_2D | H(Div) | 3 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
168 | RT0_3D | H(Div) | 1 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
169 | RT1_3D | H(Div) | 2 | 1 / 0 | H_DIV | Left in for backward compatibility, consider using RT_ |
170
171 | Tag | Description |
172 | :------: | :--------: |
173 | [DIM] | Dimension of the elements (1D, 2D, 3D) |
174 | [ORDER] | Approximation order of the elements (P0, P1, P2, ...) |
175 | [BTYPE] | BasisType of the element (0-GaussLegendre, 1 - GaussLobatto, 2-Bernstein, 3-OpenUniform, 4-CloseUniform, 5-OpenHalfUniform) |
176 | [OBTYPE] | Open BasisType of the element for elements which have both types |
177 | [CBTYPE] | Closed BasisType of the element for elements which have both types |
178
179 [FENAME] Is a special case for the Local FEC which generates a local version of a given
180 FEC. It is selected from one of (BiCubic2DFiniteElement, Quad_Q3, Nedelec1HexFiniteElement,
181 Hex_ND1, H1_[DIM]_[ORDER],H1Pos_[DIM]_[ORDER], L2_[DIM]_[ORDER] )
182 */
183 static FiniteElementCollection *New(const char *name);
184
185 /** @brief Get the local dofs for a given sub-manifold.
186
187 Return the local dofs for a SDim-dimensional sub-manifold (0D - vertex, 1D
188 - edge, 2D - face) including those on its boundary. The local index of the
189 sub-manifold (inside Geom) and its orientation are given by the parameter
190 Info = 64 * SubIndex + SubOrientation. Naturally, it is assumed that 0 <=
191 SDim <= Dim(Geom). */
192 void SubDofOrder(Geometry::Type Geom, int SDim, int Info,
193 Array<int> &dofs) const;
194
195 /// Variable order version of FiniteElementForGeometry().
196 /** The order parameter @a p represents the order of the highest-dimensional
197 FiniteElement%s the fixed-order collection we want to query. In general,
198 this order is different from the order of the returned FiniteElement. */
199 const FiniteElement *GetFE(Geometry::Type geom, int p) const
200 {
201 if (p == base_p) { return FiniteElementForGeometry(geom); }
202 if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
203 return var_orders[p]->FiniteElementForGeometry(geom);
204 }
205
206 /// Variable order version of TraceFiniteElementForGeometry().
207 /** The order parameter @a p represents the order of the highest-dimensional
208 FiniteElement%s the fixed-order collection we want to query. In general,
209 this order is different from the order of the returned FiniteElement. */
210 const FiniteElement *GetTraceFE(Geometry::Type geom, int p) const
211 {
212 if (p == base_p) { return TraceFiniteElementForGeometry(geom); }
213 if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
214 return var_orders[p]->TraceFiniteElementForGeometry(geom);
215 }
216
217 /// Variable order version of DofForGeometry().
218 /** The order parameter @a p represents the order of the highest-dimensional
219 FiniteElement%s the fixed-order collection we want to query. In general,
220 this order is different from the order of the element corresponding to
221 @a geom in that fixed-order collection. */
222 int GetNumDof(Geometry::Type geom, int p) const
223 {
224 if (p == base_p) { return DofForGeometry(geom); }
225 if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
226 return var_orders[p]->DofForGeometry(geom);
227 }
228
229 /// Variable order version of DofOrderForOrientation().
230 /** The order parameter @a p represents the order of the highest-dimensional
231 FiniteElement%s the fixed-order collection we want to query. In general,
232 this order is different from the order of the element corresponding to
233 @a geom in that fixed-order collection. */
234 const int *GetDofOrdering(Geometry::Type geom, int p, int ori) const
235 {
236 if (p == base_p) { return DofOrderForOrientation(geom, ori); }
237 if (p >= var_orders.Size() || !var_orders[p]) { InitVarOrder(p); }
238 return var_orders[p]->DofOrderForOrientation(geom, ori);
239 }
240
241 /** @brief Return the order (polynomial degree) of the FE collection,
242 corresponding to the order/degree returned by FiniteElement::GetOrder()
243 of the highest-dimensional FiniteElement%s defined by the collection. */
244 int GetOrder() const { return base_p; }
245
246 /// Instantiate a new collection of the same type with a different order.
247 /** Generally, the order parameter @a p is NOT the same as the parameter @a p
248 used by some of the constructors of derived classes. Instead, this @a p
249 represents the order of the new FE collection as it will be returned by
250 its GetOrder() method. */
251 virtual FiniteElementCollection *Clone(int p) const;
252
253protected:
254 const int base_p; ///< Order as returned by GetOrder().
255
258
259 void InitVarOrder(int p) const;
260
262
263 /// How to treat errors in FiniteElementForGeometry() calls.
265 {
266 RETURN_NULL, ///< Return NULL on errors
267 RAISE_MFEM_ERROR /**< Raise an MFEM error (default in base class).
268 Sub-classes can ignore this and return NULL. */
269 };
270
271 /// How to treat errors in FiniteElementForGeometry() calls.
272 /** The typical error in derived classes is that no FiniteElement is defined
273 for the given Geometry, or the input is not a valid Geometry. */
275};
276
277/// Arbitrary order H1-conforming (continuous) finite elements.
279{
280protected:
282 char h1_name[32];
285 int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8], *TetDofOrd[24];
286
287public:
288 explicit H1_FECollection(const int p, const int dim = 3,
289 const int btype = BasisType::GaussLobatto,
290 const int pyrtype = 1);
291
292 const FiniteElement *
293 FiniteElementForGeometry(Geometry::Type GeomType) const override;
294
295 int DofForGeometry(Geometry::Type GeomType) const override
296 { return H1_dof[GeomType]; }
297
298 const int *DofOrderForOrientation(Geometry::Type GeomType,
299 int Or) const override;
300
301 const char *Name() const override { return h1_name; }
302
303 int GetContType() const override { return CONTINUOUS; }
304
305 int GetBasisType() const { return b_type; }
306
308
309 /// Get the Cartesian to local H1 dof map
310 const int *GetDofMap(Geometry::Type GeomType) const;
311 /// Variable order version of GetDofMap
312 const int *GetDofMap(Geometry::Type GeomType, int p) const;
313
314 FiniteElementCollection *Clone(int p) const override
315 { return new H1_FECollection(p, dim, b_type); }
316
317 virtual ~H1_FECollection();
318};
319
320/** @brief Arbitrary order H1-conforming (continuous) finite elements with
321 positive basis functions. */
323{
324public:
325 explicit H1Pos_FECollection(const int p, const int dim = 3)
326 : H1_FECollection(p, dim, BasisType::Positive) {}
327};
328
329/** Arbitrary order H1-conforming (continuous) serendipity finite elements;
330 Current implementation works in 2D only; 3D version is in development. */
332{
333public:
334 explicit H1Ser_FECollection(const int p, const int dim = 2)
335 : H1_FECollection(p, dim, BasisType::Serendipity) {}
336};
337
338/** @brief Arbitrary order "H^{1/2}-conforming" trace finite elements defined on
339 the interface between mesh elements (faces,edges,vertices); these are the
340 trace FEs of the H1-conforming FEs. */
342{
343public:
344 H1_Trace_FECollection(const int p, const int dim,
345 const int btype = BasisType::GaussLobatto);
346};
347
348/// Arbitrary order "L2-conforming" discontinuous finite elements.
350{
351private:
352 int dim;
353 int b_type; // BasisType
354 int m_type; // map type
355 char d_name[32];
358 int *SegDofOrd[2]; // for rotating segment dofs in 1D
359 int *TriDofOrd[6]; // for rotating triangle dofs in 2D
360 int *TetDofOrd[24]; // for rotating tetrahedron dofs in 3D
361 int *OtherDofOrd; // for rotating other types of elements (for Or == 0)
362
363public:
364 L2_FECollection(const int p, const int dim,
365 const int btype = BasisType::GaussLegendre,
366 const int map_type = FiniteElement::VALUE,
367 const int pyrtype = 1);
368
369 const FiniteElement *
370 FiniteElementForGeometry(Geometry::Type GeomType) const override;
371
372 int DofForGeometry(Geometry::Type GeomType) const override
373 {
374 if (L2_Elements[GeomType])
375 {
376 return L2_Elements[GeomType]->GetDof();
377 }
378 return 0;
379 }
380
381 const int *DofOrderForOrientation(Geometry::Type GeomType,
382 int Or) const override;
383
384 const char *Name() const override { return d_name; }
385
386 int GetContType() const override { return DISCONTINUOUS; }
387
388 const FiniteElement *
390 {
391 return Tr_Elements[GeomType];
392 }
393
394 int GetBasisType() const { return b_type; }
395
396 FiniteElementCollection *Clone(int p) const override
397 { return new L2_FECollection(p, dim, b_type, m_type); }
398
399 virtual ~L2_FECollection();
400};
401
402/// Declare an alternative name for L2_FECollection = DG_FECollection
404
405/// Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
407{
408protected:
409 int dim;
410 int cb_type; // closed BasisType
411 int ob_type; // open BasisType
412 char rt_name[32];
415 int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
416
417 // Initialize only the face elements
418 void InitFaces(const int p, const int dim, const int map_type,
419 const bool signs);
420
421 // Constructor used by the constructor of the RT_Trace_FECollection and
422 // DG_Interface_FECollection classes
423 RT_FECollection(const int p, const int dim, const int map_type,
424 const bool signs,
426
427public:
428 /// Construct an H(div)-conforming Raviart-Thomas FE collection, RT_p.
429 /** The index @a p corresponds to the space RT_p, as typically denoted in the
430 literature, which contains vector polynomials of degree up to (p+1).
431 For example, the RT_0 collection contains vector-valued linear functions
432 and, in particular, FiniteElementCollection::GetOrder() will,
433 correspondingly, return order 1. */
434 RT_FECollection(const int p, const int dim,
437
438 const FiniteElement *
439 FiniteElementForGeometry(Geometry::Type GeomType) const override;
440
441 int DofForGeometry(Geometry::Type GeomType) const override
442 { return RT_dof[GeomType]; }
443
444 const int *DofOrderForOrientation(Geometry::Type GeomType,
445 int Or) const override;
446
447 const char *Name() const override { return rt_name; }
448
449 int GetContType() const override { return NORMAL; }
450
452
453 int GetClosedBasisType() const { return cb_type; }
454 int GetOpenBasisType() const { return ob_type; }
455
456 FiniteElementCollection *Clone(int p) const override
457 { return new RT_FECollection(p, dim, cb_type, ob_type); }
458
459 virtual ~RT_FECollection();
460};
461
462/** @brief Arbitrary order "H^{-1/2}-conforming" face finite elements defined on
463 the interface between mesh elements (faces); these are the normal trace FEs
464 of the H(div)-conforming FEs. */
466{
467public:
468 RT_Trace_FECollection(const int p, const int dim,
469 const int map_type = FiniteElement::INTEGRAL,
471
472 FiniteElementCollection *Clone(int p) const override
473 {
474 const int map_type = (strncmp(rt_name, "RT_Trace", 8) == 0)?
476 return new RT_Trace_FECollection(p, dim, map_type, ob_type);
477 }
478};
479
480/** Arbitrary order discontinuous finite elements defined on the interface
481 between mesh elements (faces). The functions in this space are single-valued
482 on each face and are discontinuous across its boundary. */
484{
485public:
486 DG_Interface_FECollection(const int p, const int dim,
487 const int map_type = FiniteElement::VALUE,
489
490 FiniteElementCollection *Clone(int p) const override
491 {
492 const int map_type = (strncmp(rt_name, "DG_Iface", 8) == 0)?
494 return new DG_Interface_FECollection(p, dim, map_type, ob_type);
495 }
496};
497
498/// Arbitrary order H(curl)-conforming Nedelec finite elements.
500{
501protected:
502 int dim;
503 int cb_type; // closed BasisType
504 int ob_type; // open BasisType
505 char nd_name[32];
508 int *SegDofOrd[2], *TriDofOrd[6], *QuadDofOrd[8];
509
510public:
511 ND_FECollection(const int p, const int dim,
514
515 const FiniteElement *
516 FiniteElementForGeometry(Geometry::Type GeomType) const override;
517
518 int DofForGeometry(Geometry::Type GeomType) const override
519 { return ND_dof[GeomType]; }
520
522 DofTransformationForGeometry(Geometry::Type GeomType) const override;
523
524 const int *DofOrderForOrientation(Geometry::Type GeomType,
525 int Or) const override;
526
527 const char *Name() const override { return nd_name; }
528
529 int GetContType() const override { return TANGENTIAL; }
530
532
533 int GetClosedBasisType() const { return cb_type; }
534 int GetOpenBasisType() const { return ob_type; }
535
536 FiniteElementCollection *Clone(int p) const override
537 { return new ND_FECollection(p, dim, cb_type, ob_type); }
538
539 virtual ~ND_FECollection();
540};
541
542/** @brief Arbitrary order H(curl)-trace finite elements defined on the
543 interface between mesh elements (faces,edges); these are the tangential
544 trace FEs of the H(curl)-conforming FEs. */
546{
547public:
548 ND_Trace_FECollection(const int p, const int dim,
551};
552
553/// Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
555{
556protected:
557 char nd_name[32];
560
561public:
562 ND_R1D_FECollection(const int p, const int dim,
563 const int cb_type = BasisType::GaussLobatto,
564 const int ob_type = BasisType::GaussLegendre);
565
566 const FiniteElement *
568 { return ND_Elements[GeomType]; }
569
570 int DofForGeometry(Geometry::Type GeomType) const override
571 { return ND_dof[GeomType]; }
572
573 const int *DofOrderForOrientation(Geometry::Type GeomType,
574 int Or) const override;
575
576 const char *Name() const override { return nd_name; }
577
578 int GetContType() const override { return TANGENTIAL; }
579
581
582 virtual ~ND_R1D_FECollection();
583};
584
585/// Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 1D.
587{
588protected:
589 char rt_name[32];
592
593public:
594 RT_R1D_FECollection(const int p, const int dim,
595 const int cb_type = BasisType::GaussLobatto,
596 const int ob_type = BasisType::GaussLegendre);
597
598 const FiniteElement *
600 { return RT_Elements[GeomType]; }
601
602 int DofForGeometry(Geometry::Type GeomType) const override
603 { return RT_dof[GeomType]; }
604
605 const int *DofOrderForOrientation(Geometry::Type GeomType,
606 int Or) const override;
607
608 const char *Name() const override { return rt_name; }
609
610 int GetContType() const override { return NORMAL; }
611
613
614 virtual ~RT_R1D_FECollection();
615};
616
617/// Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
619{
620protected:
621 char nd_name[32];
624 int *SegDofOrd[2];
625
626public:
627 ND_R2D_FECollection(const int p, const int dim,
628 const int cb_type = BasisType::GaussLobatto,
629 const int ob_type = BasisType::GaussLegendre);
630
631 const FiniteElement *
633 { return ND_Elements[GeomType]; }
634
635 int DofForGeometry(Geometry::Type GeomType) const override
636 { return ND_dof[GeomType]; }
637
638 const int *DofOrderForOrientation(Geometry::Type GeomType,
639 int Or) const override;
640
641 const char *Name() const override { return nd_name; }
642
643 int GetContType() const override { return TANGENTIAL; }
644
646
647 virtual ~ND_R2D_FECollection();
648};
649
650/** @brief Arbitrary order 3D H(curl)-trace finite elements in 2D defined on the
651 interface between mesh elements (edges); these are the tangential
652 trace FEs of the H(curl)-conforming FEs. */
654{
655public:
656 ND_R2D_Trace_FECollection(const int p, const int dim,
657 const int cb_type = BasisType::GaussLobatto,
658 const int ob_type = BasisType::GaussLegendre);
659};
660
661/// Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 2D.
663{
664protected:
665 int ob_type; // open BasisType
666 char rt_name[32];
669 int *SegDofOrd[2];
670
671 // Initialize only the face elements
672 void InitFaces(const int p, const int dim, const int map_type,
673 const bool signs);
674
675 // Constructor used by the constructor of the RT_R2D_Trace_FECollection
676 RT_R2D_FECollection(const int p, const int dim, const int map_type,
677 const bool signs,
679
680public:
681 RT_R2D_FECollection(const int p, const int dim,
682 const int cb_type = BasisType::GaussLobatto,
684
685 const FiniteElement *
687 { return RT_Elements[GeomType]; }
688
689 int DofForGeometry(Geometry::Type GeomType) const override
690 { return RT_dof[GeomType]; }
691
692 const int *DofOrderForOrientation(Geometry::Type GeomType,
693 int Or) const override;
694
695 const char *Name() const override { return rt_name; }
696
697 int GetContType() const override { return NORMAL; }
698
700
701 virtual ~RT_R2D_FECollection();
702};
703
704/** @brief Arbitrary order 3D "H^{-1/2}-conforming" face finite elements defined on
705 the interface between mesh elements (faces); these are the normal trace FEs
706 of the H(div)-conforming FEs. */
708{
709public:
710 RT_R2D_Trace_FECollection(const int p, const int dim,
711 const int map_type = FiniteElement::INTEGRAL,
713};
714
715/// Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
717{
718protected:
723
724 mutable int mOrder; // >= 1 or VariableOrder
725 // The 'name' can be:
726 // 1) name = "NURBS" + "number", for fixed order, or
727 // 2) name = "NURBS", for VariableOrder.
728 // The name is updated before writing it to a stream, for example, see
729 // FiniteElementSpace::Save().
730 mutable char name[16];
731
732public:
733 enum { VariableOrder = -1 };
734
735 /** @brief The parameter @a Order must be either a positive number, for fixed
736 order, or VariableOrder (default). */
737 explicit NURBSFECollection(int Order = VariableOrder);
738
739 virtual void Reset() const
740 {
741 SegmentFE->Reset();
744 }
745
746 virtual void SetDim(const int dim) {};
747
748 /** @brief Get the order of the NURBS collection: either a positive number,
749 when using fixed order, or VariableOrder. */
750 /** @note Not to be confused with FiniteElementCollection::GetOrder(). */
751 int GetOrder() const { return mOrder; }
752
753 /** @brief Set the order and the name, based on the given @a Order: either a
754 positive number for fixed order, or VariableOrder. */
755 virtual void SetOrder(int Order) const;
756
757 const FiniteElement *
758 FiniteElementForGeometry(Geometry::Type GeomType) const override;
759
760 int DofForGeometry(Geometry::Type GeomType) const override;
761
762 const int *DofOrderForOrientation(Geometry::Type GeomType,
763 int Or) const override;
764
765 const char *Name() const override { return name; }
766
767 int GetContType() const override { return CONTINUOUS; }
768
770
771 virtual ~NURBSFECollection();
772};
773
774/// Arbitrary order H(div) NURBS finite elements.
776{
777private:
778
779 NURBS1DFiniteElement *SegmentFE;
780 NURBS2DFiniteElement *QuadrilateralFE;
781
782 NURBS_HDiv2DFiniteElement *QuadrilateralVFE;
783 NURBS_HDiv3DFiniteElement *ParallelepipedVFE;
784
785 FiniteElement *sFE;
786 FiniteElement *qFE;
787 FiniteElement *hFE;
788
789public:
790
791 /** @brief The parameter @a Order must be either a positive number, for fixed
792 order, or VariableOrder (default). */
793 explicit NURBS_HDivFECollection(int Order = VariableOrder, const int vdim = -1);
794
795 void Reset() const override
796 {
797 SegmentFE->Reset();
798 QuadrilateralFE->Reset();
799 QuadrilateralVFE->Reset();
800 ParallelepipedVFE->Reset();
801 }
802
803 void SetDim(const int dim) override;
804
805 /** @brief Set the order and the name, based on the given @a Order: either a
806 positive number for fixed order, or VariableOrder. */
807 void SetOrder(int Order) const override;
808
809 const FiniteElement *
810 FiniteElementForGeometry(Geometry::Type GeomType) const override;
811
812 int DofForGeometry(Geometry::Type GeomType) const override;
813
814 const int *DofOrderForOrientation(Geometry::Type GeomType,
815 int Or) const override;
816
817 const char *Name() const override { return name; }
818
819 int GetContType() const override { return CONTINUOUS; }
820
822
823 virtual ~NURBS_HDivFECollection();
824};
825
826/// Arbitrary order H(curl) NURBS finite elements.
828{
829private:
830 NURBS1DFiniteElement *SegmentFE;
831 NURBS2DFiniteElement *QuadrilateralFE;
832
833 NURBS_HCurl2DFiniteElement *QuadrilateralVFE;
834 NURBS_HCurl3DFiniteElement *ParallelepipedVFE;
835
836 FiniteElement *sFE;
837 FiniteElement *qFE;
838 FiniteElement *hFE;
839public:
840
841 /** @brief The parameter @a Order must be either a positive number, for fixed
842 order, or VariableOrder (default). */
843 explicit NURBS_HCurlFECollection(int Order = VariableOrder,
844 const int vdim = -1);
845
846 void Reset() const override
847 {
848 SegmentFE->Reset();
849 QuadrilateralFE->Reset();
850 QuadrilateralVFE->Reset();
851 ParallelepipedVFE->Reset();
852 }
853
854 void SetDim(const int dim) override;
855
856 /** @brief Set the order and the name, based on the given @a Order: either a
857 positive number for fixed order, or VariableOrder. */
858 void SetOrder(int Order) const override;
859
860 const FiniteElement *
861 FiniteElementForGeometry(Geometry::Type GeomType) const override;
862
863 int DofForGeometry(Geometry::Type GeomType) const override;
864
865 const int *DofOrderForOrientation(Geometry::Type GeomType,
866 int Or) const override;
867
868 const char *Name() const override { return name; }
869
870 int GetContType() const override { return CONTINUOUS; }
871
873
874 virtual ~NURBS_HCurlFECollection();
875};
876
877/// Piecewise-(bi/tri)linear continuous finite elements.
879{
880private:
881 const PointFiniteElement PointFE;
882 const Linear1DFiniteElement SegmentFE;
883 const Linear2DFiniteElement TriangleFE;
884 const BiLinear2DFiniteElement QuadrilateralFE;
885 const Linear3DFiniteElement TetrahedronFE;
886 const TriLinear3DFiniteElement ParallelepipedFE;
887 const LinearWedgeFiniteElement WedgeFE;
888 const LinearPyramidFiniteElement PyramidFE;
889public:
891
892 const FiniteElement *
893 FiniteElementForGeometry(Geometry::Type GeomType) const override;
894
895 int DofForGeometry(Geometry::Type GeomType) const override;
896
897 const int *DofOrderForOrientation(Geometry::Type GeomType,
898 int Or) const override;
899
900 const char *Name() const override { return "Linear"; }
901
902 int GetContType() const override { return CONTINUOUS; }
903};
904
905/// Piecewise-(bi)quadratic continuous finite elements.
907{
908private:
909 const PointFiniteElement PointFE;
910 const Quad1DFiniteElement SegmentFE;
911 const Quad2DFiniteElement TriangleFE;
912 const BiQuad2DFiniteElement QuadrilateralFE;
913 const Quadratic3DFiniteElement TetrahedronFE;
914 const LagrangeHexFiniteElement ParallelepipedFE;
915 const H1_WedgeElement WedgeFE;
916
917public:
919 : FiniteElementCollection(2), ParallelepipedFE(2), WedgeFE(2) {}
920
921 const FiniteElement *
922 FiniteElementForGeometry(Geometry::Type GeomType) const override;
923
924 int DofForGeometry(Geometry::Type GeomType) const override;
925
926 const int *DofOrderForOrientation(Geometry::Type GeomType,
927 int Or) const override;
928
929 const char *Name() const override { return "Quadratic"; }
930
931 int GetContType() const override { return CONTINUOUS; }
932};
933
934/// Version of QuadraticFECollection with positive basis functions.
936{
937private:
938 const QuadPos1DFiniteElement SegmentFE;
939 const BiQuadPos2DFiniteElement QuadrilateralFE;
940
941public:
943
944 const FiniteElement *
945 FiniteElementForGeometry(Geometry::Type GeomType) const override;
946
947 int DofForGeometry(Geometry::Type GeomType) const override;
948
949 const int *DofOrderForOrientation(Geometry::Type GeomType,
950 int Or) const override;
951
952 const char *Name() const override { return "QuadraticPos"; }
953
954 int GetContType() const override { return CONTINUOUS; }
955};
956
957/// Piecewise-(bi)cubic continuous finite elements.
959{
960private:
961 const PointFiniteElement PointFE;
962 const Cubic1DFiniteElement SegmentFE;
963 const Cubic2DFiniteElement TriangleFE;
964 const BiCubic2DFiniteElement QuadrilateralFE;
965 const Cubic3DFiniteElement TetrahedronFE;
966 const LagrangeHexFiniteElement ParallelepipedFE;
967 const H1_WedgeElement WedgeFE;
968
969public:
972 ParallelepipedFE(3), WedgeFE(3, BasisType::ClosedUniform)
973 {}
974
975 const FiniteElement *
976 FiniteElementForGeometry(Geometry::Type GeomType) const override;
977
978 int DofForGeometry(Geometry::Type GeomType) const override;
979
980 const int *DofOrderForOrientation(Geometry::Type GeomType,
981 int Or) const override;
982
983 const char *Name() const override { return "Cubic"; }
984
985 int GetContType() const override { return CONTINUOUS; }
986};
987
988/// Crouzeix-Raviart nonconforming elements in 2D.
990{
991private:
992 const P0SegmentFiniteElement SegmentFE;
993 const CrouzeixRaviartFiniteElement TriangleFE;
994 const CrouzeixRaviartQuadFiniteElement QuadrilateralFE;
995public:
997
998 const FiniteElement *
999 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1000
1001 int DofForGeometry(Geometry::Type GeomType) const override;
1002
1003 const int *DofOrderForOrientation(Geometry::Type GeomType,
1004 int Or) const override;
1005
1006 const char *Name() const override { return "CrouzeixRaviart"; }
1007
1008 int GetContType() const override { return DISCONTINUOUS; }
1009};
1010
1011/// Piecewise-linear nonconforming finite elements in 3D.
1013{
1014private:
1015 const P0TriangleFiniteElement TriangleFE;
1016 const P1TetNonConfFiniteElement TetrahedronFE;
1017 const P0QuadFiniteElement QuadrilateralFE;
1018 const RotTriLinearHexFiniteElement ParallelepipedFE;
1019
1020public:
1022
1023 const FiniteElement *
1024 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1025
1026 int DofForGeometry(Geometry::Type GeomType) const override;
1027
1028 const int *DofOrderForOrientation(Geometry::Type GeomType,
1029 int Or) const override;
1030
1031 const char *Name() const override { return "LinearNonConf3D"; }
1032
1033 int GetContType() const override { return DISCONTINUOUS; }
1034};
1035
1036/** @brief First order Raviart-Thomas finite elements in 2D. This class is kept
1037 only for backward compatibility, consider using RT_FECollection instead. */
1039{
1040private:
1041 const P0SegmentFiniteElement SegmentFE; // normal component on edge
1042 const RT0TriangleFiniteElement TriangleFE;
1043 const RT0QuadFiniteElement QuadrilateralFE;
1044public:
1046
1047 const FiniteElement *
1048 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1049
1050 int DofForGeometry(Geometry::Type GeomType) const override;
1051
1052 const int *DofOrderForOrientation(Geometry::Type GeomType,
1053 int Or) const override;
1054
1055 const char *Name() const override { return "RT0_2D"; }
1056
1057 int GetContType() const override { return NORMAL; }
1058};
1059
1060/** @brief Second order Raviart-Thomas finite elements in 2D. This class is kept
1061 only for backward compatibility, consider using RT_FECollection instead. */
1063{
1064private:
1065 const P1SegmentFiniteElement SegmentFE; // normal component on edge
1066 const RT1TriangleFiniteElement TriangleFE;
1067 const RT1QuadFiniteElement QuadrilateralFE;
1068public:
1070
1071 const FiniteElement *
1072 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1073
1074 int DofForGeometry(Geometry::Type GeomType) const override;
1075
1076 const int *DofOrderForOrientation(Geometry::Type GeomType,
1077 int Or) const override;
1078
1079 const char *Name() const override { return "RT1_2D"; }
1080
1081 int GetContType() const override { return NORMAL; }
1082};
1083
1084/** @brief Third order Raviart-Thomas finite elements in 2D. This class is kept
1085 only for backward compatibility, consider using RT_FECollection instead. */
1087{
1088private:
1089 const P2SegmentFiniteElement SegmentFE; // normal component on edge
1090 const RT2TriangleFiniteElement TriangleFE;
1091 const RT2QuadFiniteElement QuadrilateralFE;
1092public:
1094
1095 const FiniteElement *
1096 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1097
1098 int DofForGeometry(Geometry::Type GeomType) const override;
1099
1100 const int *DofOrderForOrientation(Geometry::Type GeomType,
1101 int Or) const override;
1102
1103 const char *Name() const override { return "RT2_2D"; }
1104
1105 int GetContType() const override { return NORMAL; }
1106};
1107
1108/** @brief Piecewise-constant discontinuous finite elements in 2D. This class is
1109 kept only for backward compatibility, consider using L2_FECollection
1110 instead. */
1112{
1113private:
1114 const P0TriangleFiniteElement TriangleFE;
1115 const P0QuadFiniteElement QuadrilateralFE;
1116public:
1118
1119 const FiniteElement *
1120 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1121
1122 int DofForGeometry(Geometry::Type GeomType) const override;
1123
1124 const int *DofOrderForOrientation(Geometry::Type GeomType,
1125 int Or) const override;
1126
1127 const char *Name() const override { return "Const2D"; }
1128
1129 int GetContType() const override { return DISCONTINUOUS; }
1130};
1131
1132/** @brief Piecewise-linear discontinuous finite elements in 2D. This class is
1133 kept only for backward compatibility, consider using L2_FECollection
1134 instead. */
1136{
1137private:
1138 const Linear2DFiniteElement TriangleFE;
1139 const BiLinear2DFiniteElement QuadrilateralFE;
1140
1141public:
1143
1144 const FiniteElement *
1145 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1146
1147 int DofForGeometry(Geometry::Type GeomType) const override;
1148
1149 const int *DofOrderForOrientation(Geometry::Type GeomType,
1150 int Or) const override;
1151
1152 const char *Name() const override { return "LinearDiscont2D"; }
1153
1154 int GetContType() const override { return DISCONTINUOUS; }
1155};
1156
1157/// Version of LinearDiscont2DFECollection with dofs in the Gaussian points.
1159{
1160private:
1161 // const CrouzeixRaviartFiniteElement TriangleFE;
1162 const GaussLinear2DFiniteElement TriangleFE;
1163 const GaussBiLinear2DFiniteElement QuadrilateralFE;
1164
1165public:
1167
1168 const FiniteElement *
1169 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1170
1171 int DofForGeometry(Geometry::Type GeomType) const override;
1172
1173 const int *DofOrderForOrientation(Geometry::Type GeomType,
1174 int Or) const override;
1175
1176 const char *Name() const override { return "GaussLinearDiscont2D"; }
1177
1178 int GetContType() const override { return DISCONTINUOUS; }
1179};
1180
1181/// Linear (P1) finite elements on quadrilaterals.
1183{
1184private:
1185 const P1OnQuadFiniteElement QuadrilateralFE;
1186public:
1188
1189 const FiniteElement *
1190 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1191
1192 int DofForGeometry(Geometry::Type GeomType) const override;
1193
1194 const int *DofOrderForOrientation(Geometry::Type GeomType,
1195 int Or) const override;
1196
1197 const char *Name() const override { return "P1OnQuad"; }
1198
1199 int GetContType() const override { return DISCONTINUOUS; }
1200};
1201
1202/** @brief Piecewise-quadratic discontinuous finite elements in 2D. This class
1203 is kept only for backward compatibility, consider using L2_FECollection
1204 instead. */
1206{
1207private:
1208 const Quad2DFiniteElement TriangleFE;
1209 const BiQuad2DFiniteElement QuadrilateralFE;
1210
1211public:
1213
1214 const FiniteElement *
1215 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1216
1217 int DofForGeometry(Geometry::Type GeomType) const override;
1218
1219 const int *DofOrderForOrientation(Geometry::Type GeomType,
1220 int Or) const override;
1221
1222 const char *Name() const override { return "QuadraticDiscont2D"; }
1223
1224 int GetContType() const override { return DISCONTINUOUS; }
1225};
1226
1227/// Version of QuadraticDiscont2DFECollection with positive basis functions.
1229{
1230private:
1231 const BiQuadPos2DFiniteElement QuadrilateralFE;
1232
1233public:
1235
1236 const FiniteElement *
1237 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1238
1239 int DofForGeometry(Geometry::Type GeomType) const override;
1240
1242 int Or) const override
1243 { return NULL; }
1244
1245 const char *Name() const override { return "QuadraticPosDiscont2D"; }
1246
1247 int GetContType() const override { return DISCONTINUOUS; }
1248};
1249
1250/// Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
1252{
1253private:
1254 // const Quad2DFiniteElement TriangleFE;
1255 const GaussQuad2DFiniteElement TriangleFE;
1256 const GaussBiQuad2DFiniteElement QuadrilateralFE;
1257
1258public:
1260
1261 const FiniteElement *
1262 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1263
1264 int DofForGeometry(Geometry::Type GeomType) const override;
1265
1266 const int *DofOrderForOrientation(Geometry::Type GeomType,
1267 int Or) const override;
1268
1269 const char *Name() const override { return "GaussQuadraticDiscont2D"; }
1270
1271 int GetContType() const override { return DISCONTINUOUS; }
1272};
1273
1274/** @brief Piecewise-cubic discontinuous finite elements in 2D. This class is
1275 kept only for backward compatibility, consider using L2_FECollection
1276 instead. */
1278{
1279private:
1280 const Cubic2DFiniteElement TriangleFE;
1281 const BiCubic2DFiniteElement QuadrilateralFE;
1282
1283public:
1285
1286 const FiniteElement *
1287 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1288
1289 int DofForGeometry(Geometry::Type GeomType) const override;
1290
1291 const int *DofOrderForOrientation(Geometry::Type GeomType,
1292 int Or) const override;
1293
1294 const char *Name() const override { return "CubicDiscont2D"; }
1295
1296 int GetContType() const override { return DISCONTINUOUS; }
1297};
1298
1299/** @brief Piecewise-constant discontinuous finite elements in 3D. This class is
1300 kept only for backward compatibility, consider using L2_FECollection
1301 instead. */
1303{
1304private:
1305 const P0TetFiniteElement TetrahedronFE;
1306 const P0HexFiniteElement ParallelepipedFE;
1307 const P0WdgFiniteElement WedgeFE;
1308 const P0PyrFiniteElement PyramidFE;
1309
1310public:
1312
1313 const FiniteElement *
1314 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1315
1316 int DofForGeometry(Geometry::Type GeomType) const override;
1317
1318 const int *DofOrderForOrientation(Geometry::Type GeomType,
1319 int Or) const override;
1320
1321 const char *Name() const override { return "Const3D"; }
1322
1323 int GetContType() const override { return DISCONTINUOUS; }
1324};
1325
1326/** @brief Piecewise-linear discontinuous finite elements in 3D. This class is
1327 kept only for backward compatibility, consider using L2_FECollection
1328 instead. */
1330{
1331private:
1332 const Linear3DFiniteElement TetrahedronFE;
1333 const LinearPyramidFiniteElement PyramidFE;
1334 const LinearWedgeFiniteElement WedgeFE;
1335 const TriLinear3DFiniteElement ParallelepipedFE;
1336
1337public:
1339
1340 const FiniteElement *
1341 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1342
1343 int DofForGeometry(Geometry::Type GeomType) const override;
1344
1345 const int *DofOrderForOrientation(Geometry::Type GeomType,
1346 int Or) const override;
1347
1348 const char *Name() const override { return "LinearDiscont3D"; }
1349
1350 int GetContType() const override { return DISCONTINUOUS; }
1351};
1352
1353/** @brief Piecewise-quadratic discontinuous finite elements in 3D. This class
1354 is kept only for backward compatibility, consider using L2_FECollection
1355 instead. */
1357{
1358private:
1359 const Quadratic3DFiniteElement TetrahedronFE;
1360 const LagrangeHexFiniteElement ParallelepipedFE;
1361
1362public:
1364 : FiniteElementCollection(2), ParallelepipedFE(2) {}
1365
1366 const FiniteElement *
1367 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1368
1369 int DofForGeometry(Geometry::Type GeomType) const override;
1370
1371 const int *DofOrderForOrientation(Geometry::Type GeomType,
1372 int Or) const override;
1373
1374 const char *Name() const override { return "QuadraticDiscont3D"; }
1375
1376 int GetContType() const override { return DISCONTINUOUS; }
1377};
1378
1379/// Finite element collection on a macro-element.
1381{
1382private:
1383 const PointFiniteElement PointFE;
1384 const RefinedLinear1DFiniteElement SegmentFE;
1385 const RefinedLinear2DFiniteElement TriangleFE;
1386 const RefinedBiLinear2DFiniteElement QuadrilateralFE;
1387 const RefinedLinear3DFiniteElement TetrahedronFE;
1388 const RefinedTriLinear3DFiniteElement ParallelepipedFE;
1389
1390public:
1392
1393 const FiniteElement *
1394 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1395
1396 int DofForGeometry(Geometry::Type GeomType) const override;
1397
1398 const int *DofOrderForOrientation(Geometry::Type GeomType,
1399 int Or) const override;
1400
1401 const char *Name() const override { return "RefinedLinear"; }
1402
1403 int GetContType() const override { return CONTINUOUS; }
1404};
1405
1406/** @brief Lowest order Nedelec finite elements in 3D. This class is kept only
1407 for backward compatibility, consider using the new ND_FECollection
1408 instead. */
1410{
1411private:
1412 const Nedelec1HexFiniteElement HexahedronFE;
1413 const Nedelec1TetFiniteElement TetrahedronFE;
1414 const Nedelec1WdgFiniteElement WedgeFE;
1415 const Nedelec1PyrFiniteElement PyramidFE;
1416
1417public:
1419
1420 const FiniteElement *
1421 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1422
1423 int DofForGeometry(Geometry::Type GeomType) const override;
1424
1425 const int *DofOrderForOrientation(Geometry::Type GeomType,
1426 int Or) const override;
1427
1428 const char *Name() const override { return "ND1_3D"; }
1429
1430 int GetContType() const override { return TANGENTIAL; }
1431};
1432
1433/** @brief First order Raviart-Thomas finite elements in 3D. This class is kept
1434 only for backward compatibility, consider using RT_FECollection instead. */
1436{
1437private:
1438 const P0TriangleFiniteElement TriangleFE;
1439 const P0QuadFiniteElement QuadrilateralFE;
1440 const RT0HexFiniteElement HexahedronFE;
1441 const RT0TetFiniteElement TetrahedronFE;
1442 const RT0WdgFiniteElement WedgeFE;
1443 const RT0PyrFiniteElement PyramidFE;
1444public:
1446
1447 const FiniteElement *
1448 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1449
1450 int DofForGeometry(Geometry::Type GeomType) const override;
1451
1452 const int *DofOrderForOrientation(Geometry::Type GeomType,
1453 int Or) const override;
1454
1455 const char *Name() const override { return "RT0_3D"; }
1456
1457 int GetContType() const override { return NORMAL; }
1458};
1459
1460/** @brief Second order Raviart-Thomas finite elements in 3D. This class is kept
1461 only for backward compatibility, consider using RT_FECollection instead. */
1463{
1464private:
1465 const Linear2DFiniteElement TriangleFE;
1466 const BiLinear2DFiniteElement QuadrilateralFE;
1467 const RT1HexFiniteElement HexahedronFE;
1468public:
1470
1471 const FiniteElement *
1472 FiniteElementForGeometry(Geometry::Type GeomType) const override;
1473
1474 int DofForGeometry(Geometry::Type GeomType) const override;
1475
1476 const int *DofOrderForOrientation(Geometry::Type GeomType,
1477 int Or) const override;
1478
1479 const char *Name() const override { return "RT1_3D"; }
1480
1481 int GetContType() const override { return NORMAL; }
1482};
1483
1484/// Discontinuous collection defined locally by a given finite element.
1486{
1487private:
1488 char d_name[32];
1489 Geometry::Type GeomType;
1490 FiniteElement *Local_Element;
1491
1492public:
1493 Local_FECollection(const char *fe_name);
1494
1495 const FiniteElement *
1497 { return (GeomType == GeomType_) ? Local_Element : NULL; }
1498
1499 int DofForGeometry(Geometry::Type GeomType_) const override
1500 { return (GeomType == GeomType_) ? Local_Element->GetDof() : 0; }
1501
1503 int Or) const override
1504 { return NULL; }
1505
1506 const char *Name() const override { return d_name; }
1507
1508 int GetContType() const override { return DISCONTINUOUS; }
1509
1510 virtual ~Local_FECollection() { delete Local_Element; }
1511};
1512
1513}
1514
1515#endif
Possible basis types. Note that not all elements can use all BasisType(s).
Definition fe_base.hpp:30
@ GaussLobatto
Closed type.
Definition fe_base.hpp:36
@ GaussLegendre
Open type.
Definition fe_base.hpp:35
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:1112
const char * Name() const override
Definition fe_coll.hpp:1127
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1082
int GetContType() const override
Definition fe_coll.hpp:1129
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:1096
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1069
Piecewise-constant discontinuous finite elements in 3D. This class is kept only for backward compatib...
Definition fe_coll.hpp:1303
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1409
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1394
int GetContType() const override
Definition fe_coll.hpp:1323
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:1427
const char * Name() const override
Definition fe_coll.hpp:1321
Crouzeix-Raviart nonconforming elements in 2D.
Definition fe_coll.hpp:990
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:917
const char * Name() const override
Definition fe_coll.hpp:1006
int GetContType() const override
Definition fe_coll.hpp:1008
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:931
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:902
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:1278
const char * Name() const override
Definition fe_coll.hpp:1294
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:1344
int GetContType() const override
Definition fe_coll.hpp:1296
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1330
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1316
Piecewise-(bi)cubic continuous finite elements.
Definition fe_coll.hpp:959
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:850
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:867
const char * Name() const override
Definition fe_coll.hpp:983
int GetContType() const override
Definition fe_coll.hpp:985
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:832
FiniteElementCollection * Clone(int p) const override
Instantiate a new collection of the same type with a different order.
Definition fe_coll.hpp:490
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:2745
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:274
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:124
const int * GetDofOrdering(Geometry::Type geom, int p, int ori) const
Variable order version of DofOrderForOrientation().
Definition fe_coll.hpp:234
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition fe_coll.hpp:244
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:432
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:222
virtual const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const
Definition fe_coll.hpp:97
const FiniteElement * GetTraceFE(Geometry::Type geom, int p) const
Variable order version of TraceFiniteElementForGeometry().
Definition fe_coll.hpp:210
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:468
void InitVarOrder(int p) const
Definition fe_coll.cpp:440
int GetDerivType(int dim) const
Definition fe_coll.cpp:70
const int base_p
Order as returned by GetOrder().
Definition fe_coll.hpp:254
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:488
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:199
@ 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
Array< FiniteElementCollection * > var_orders
Definition fe_coll.hpp:261
virtual const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const =0
virtual FiniteElementCollection * GetTraceCollection() const
Definition fe_coll.cpp:118
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:533
static void GetNVE(int &nv, int &ne)
Definition fe_coll.cpp:458
ErrorMode
How to treat errors in FiniteElementForGeometry() calls.
Definition fe_coll.hpp:265
@ RETURN_NULL
Return NULL on errors.
Definition fe_coll.hpp:266
Abstract class for all finite elements.
Definition fe_base.hpp:247
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:337
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:1159
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:1171
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1140
const char * Name() const override
Definition fe_coll.hpp:1176
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1155
A quadratic element on triangle with nodes at the "Gaussian" points.
Version of QuadraticDiscont2DFECollection with dofs in the Gaussian points.
Definition fe_coll.hpp:1252
const char * Name() const override
Definition fe_coll.hpp:1269
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1292
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1276
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:1308
static const int NumGeom
Definition geom.hpp:46
Arbitrary order H1-conforming (continuous) finite elements with positive basis functions.
Definition fe_coll.hpp:323
H1Pos_FECollection(const int p, const int dim=3)
Definition fe_coll.hpp:325
H1Ser_FECollection(const int p, const int dim=2)
Definition fe_coll.hpp:334
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:279
FiniteElementCollection * Clone(int p) const override
Instantiate a new collection of the same type with a different order.
Definition fe_coll.hpp:314
H1_FECollection(const int p, const int dim=3, const int btype=BasisType::GaussLobatto, const int pyrtype=1)
Definition fe_coll.cpp:1711
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:2042
virtual ~H1_FECollection()
Definition fe_coll.cpp:2113
int GetBasisType() const
Definition fe_coll.hpp:305
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:295
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:2069
const char * Name() const override
Definition fe_coll.hpp:301
int GetContType() const override
Definition fe_coll.hpp:303
int H1_dof[Geometry::NumGeom]
Definition fe_coll.hpp:284
FiniteElement * H1_Elements[Geometry::NumGeom]
Definition fe_coll.hpp:283
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:2047
const int * GetDofMap(Geometry::Type GeomType) const
Get the Cartesian to local H1 dof map.
Definition fe_coll.cpp:2088
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
Definition fe_coll.hpp:342
H1_Trace_FECollection(const int p, const int dim, const int btype=BasisType::GaussLobatto)
Definition fe_coll.cpp:2126
Arbitrary order H1 elements in 3D on a wedge.
Definition fe_h1.hpp:132
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:350
const FiniteElement * TraceFiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:389
int GetBasisType() const
Definition fe_coll.hpp:394
L2_FECollection(const int p, const int dim, const int btype=BasisType::GaussLegendre, const int map_type=FiniteElement::VALUE, const int pyrtype=1)
Definition fe_coll.cpp:2146
virtual ~L2_FECollection()
Definition fe_coll.cpp:2451
const char * Name() const override
Definition fe_coll.hpp:384
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:2427
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:372
int GetContType() const override
Definition fe_coll.hpp:386
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:2432
FiniteElementCollection * Clone(int p) const override
Instantiate a new collection of the same type with a different order.
Definition fe_coll.hpp:396
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:1136
int GetContType() const override
Definition fe_coll.hpp:1154
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1104
int DofForGeometry(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:1132
const char * Name() const override
Definition fe_coll.hpp:1152
Piecewise-linear discontinuous finite elements in 3D. This class is kept only for backward compatibil...
Definition fe_coll.hpp:1330
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:1469
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1435
const char * Name() const override
Definition fe_coll.hpp:1348
int GetContType() const override
Definition fe_coll.hpp:1350
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1451
Piecewise-(bi/tri)linear continuous finite elements.
Definition fe_coll.hpp:879
int GetContType() const override
Definition fe_coll.hpp:902
const char * Name() const override
Definition fe_coll.hpp:900
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:724
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:705
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:742
Piecewise-linear nonconforming finite elements in 3D.
Definition fe_coll.hpp:1013
int GetContType() const override
Definition fe_coll.hpp:1033
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1368
const char * Name() const override
Definition fe_coll.hpp:1031
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:1384
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1352
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:1486
int GetContType() const override
Definition fe_coll.hpp:1508
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType_) const override
Definition fe_coll.hpp:1496
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:1502
int DofForGeometry(Geometry::Type GeomType_) const override
Definition fe_coll.hpp:1499
const char * Name() const override
Definition fe_coll.hpp:1506
Local_FECollection(const char *fe_name)
Definition fe_coll.cpp:3480
Lowest order Nedelec finite elements in 3D. This class is kept only for backward compatibility,...
Definition fe_coll.hpp:1410
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1574
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1559
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:1592
int GetContType() const override
Definition fe_coll.hpp:1430
const char * Name() const override
Definition fe_coll.hpp:1428
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition fe_coll.hpp:500
int ND_dof[Geometry::NumGeom]
Definition fe_coll.hpp:507
const char * Name() const override
Definition fe_coll.hpp:527
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:2948
FiniteElementCollection * Clone(int p) const override
Instantiate a new collection of the same type with a different order.
Definition fe_coll.hpp:536
int GetOpenBasisType() const
Definition fe_coll.hpp:534
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:2966
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:2984
virtual ~ND_FECollection()
Definition fe_coll.cpp:3004
int GetClosedBasisType() const
Definition fe_coll.hpp:533
int GetContType() const override
Definition fe_coll.hpp:529
const StatelessDofTransformation * DofTransformationForGeometry(Geometry::Type GeomType) const override
Returns a DoF transformation object compatible with this basis and geometry type.
Definition fe_coll.cpp:2954
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition fe_coll.hpp:506
ND_FECollection(const int p, const int dim, const int cb_type=BasisType::GaussLobatto, const int ob_type=BasisType::GaussLegendre)
Definition fe_coll.cpp:2765
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:518
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
Definition fe_coll.hpp:555
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:3084
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:570
int ND_dof[Geometry::NumGeom]
Definition fe_coll.hpp:559
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:567
int GetContType() const override
Definition fe_coll.hpp:578
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:3035
const char * Name() const override
Definition fe_coll.hpp:576
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition fe_coll.hpp:558
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:3090
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
Definition fe_coll.hpp:619
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:3267
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:3174
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:632
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:635
FiniteElement * ND_Elements[Geometry::NumGeom]
Definition fe_coll.hpp:622
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:3257
int GetContType() const override
Definition fe_coll.hpp:643
const char * Name() const override
Definition fe_coll.hpp:641
int ND_dof[Geometry::NumGeom]
Definition fe_coll.hpp:623
Arbitrary order 3D H(curl)-trace finite elements in 2D defined on the interface between mesh elements...
Definition fe_coll.hpp:654
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:3297
Arbitrary order H(curl)-trace finite elements defined on the interface between mesh elements (faces,...
Definition fe_coll.hpp:546
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:3016
An arbitrary order 1D NURBS element on a segment.
Definition fe_nurbs.hpp:72
An arbitrary order 2D NURBS element on a square.
Definition fe_nurbs.hpp:94
An arbitrary order 3D NURBS element on a cube.
Definition fe_nurbs.hpp:129
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition fe_coll.hpp:717
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:3571
int GetOrder() const
Get the order of the NURBS collection: either a positive number, when using fixed order,...
Definition fe_coll.hpp:751
virtual void Reset() const
Definition fe_coll.hpp:739
virtual 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:3534
const char * Name() const override
Definition fe_coll.hpp:765
virtual void SetDim(const int dim)
Definition fe_coll.hpp:746
int GetContType() const override
Definition fe_coll.hpp:767
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:3556
NURBS2DFiniteElement * QuadrilateralFE
Definition fe_coll.hpp:721
NURBS1DFiniteElement * SegmentFE
Definition fe_coll.hpp:720
NURBSFECollection(int Order=VariableOrder)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default).
Definition fe_coll.cpp:3522
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:3584
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:3577
NURBS3DFiniteElement * ParallelepipedFE
Definition fe_coll.hpp:722
PointFiniteElement * PointFE
Definition fe_coll.hpp:719
void Reset() const
Resets the patch and element data stored in the element.
Definition fe_nurbs.hpp:43
Arbitrary order H(curl) NURBS finite elements.
Definition fe_coll.hpp:828
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:3762
void SetOrder(int Order) const override
Set the order and the name, based on the given Order: either a positive number for fixed order,...
Definition fe_coll.cpp:3743
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:3729
const char * Name() const override
Definition fe_coll.hpp:868
NURBS_HCurlFECollection(int Order=VariableOrder, const int vdim=-1)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default).
Definition fe_coll.cpp:3683
void SetDim(const int dim) override
Definition fe_coll.cpp:3697
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:3756
void Reset() const override
Definition fe_coll.hpp:846
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:3770
int GetContType() const override
Definition fe_coll.hpp:870
Arbitrary order H(div) NURBS finite elements.
Definition fe_coll.hpp:776
void SetDim(const int dim) override
Definition fe_coll.cpp:3606
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:3663
void SetOrder(int Order) const override
Set the order and the name, based on the given Order: either a positive number for fixed order,...
Definition fe_coll.cpp:3650
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:3677
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:3669
int GetContType() const override
Definition fe_coll.hpp:819
const char * Name() const override
Definition fe_coll.hpp:817
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:3636
void Reset() const override
Definition fe_coll.hpp:795
NURBS_HDivFECollection(int Order=VariableOrder, const int vdim=-1)
The parameter Order must be either a positive number, for fixed order, or VariableOrder (default).
Definition fe_coll.cpp:3591
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:1183
const char * Name() const override
Definition fe_coll.hpp:1197
int GetContType() const override
Definition fe_coll.hpp:1199
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1179
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:1202
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1189
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:109
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:1206
const char * Name() const override
Definition fe_coll.hpp:1222
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1224
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1210
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:1239
Piecewise-quadratic discontinuous finite elements in 3D. This class is kept only for backward compati...
Definition fe_coll.hpp:1357
const char * Name() const override
Definition fe_coll.hpp:1374
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1477
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1491
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:1508
Piecewise-(bi)quadratic continuous finite elements.
Definition fe_coll.hpp:907
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:768
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:785
int GetContType() const override
Definition fe_coll.hpp:931
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:750
const char * Name() const override
Definition fe_coll.hpp:929
Version of QuadraticDiscont2DFECollection with positive basis functions.
Definition fe_coll.hpp:1229
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1247
const char * Name() const override
Definition fe_coll.hpp:1245
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1260
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:1241
Version of QuadraticFECollection with positive basis functions.
Definition fe_coll.hpp:936
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:795
int GetContType() const override
Definition fe_coll.hpp:954
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:822
const char * Name() const override
Definition fe_coll.hpp:952
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:809
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:1039
const char * Name() const override
Definition fe_coll.hpp:1055
int GetContType() const override
Definition fe_coll.hpp:1057
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:941
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:969
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:955
First order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility,...
Definition fe_coll.hpp:1436
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:1642
const char * Name() const override
Definition fe_coll.hpp:1455
int GetContType() const override
Definition fe_coll.hpp:1457
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1607
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1624
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:1063
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:984
const char * Name() const override
Definition fe_coll.hpp:1079
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:1012
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:998
int GetContType() const override
Definition fe_coll.hpp:1081
Second order Raviart-Thomas finite elements in 3D. This class is kept only for backward compatibility...
Definition fe_coll.hpp:1463
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1674
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1660
int GetContType() const override
Definition fe_coll.hpp:1481
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:1689
const char * Name() const override
Definition fe_coll.hpp:1479
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:1087
const char * Name() const override
Definition fe_coll.hpp:1103
int GetContType() const override
Definition fe_coll.hpp:1105
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1040
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:1054
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1026
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:407
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:441
int RT_dof[Geometry::NumGeom]
Definition fe_coll.hpp:414
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition fe_coll.hpp:413
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition fe_coll.cpp:2553
int GetClosedBasisType() const
Definition fe_coll.hpp:453
virtual ~RT_FECollection()
Definition fe_coll.cpp:2713
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:2678
const char * Name() const override
Definition fe_coll.hpp:447
FiniteElementCollection * Clone(int p) const override
Instantiate a new collection of the same type with a different order.
Definition fe_coll.hpp:456
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:2673
int GetContType() const override
Definition fe_coll.hpp:449
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:2696
int GetOpenBasisType() const
Definition fe_coll.hpp:454
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:2537
Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 1D.
Definition fe_coll.hpp:587
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:3159
int GetContType() const override
Definition fe_coll.hpp:610
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition fe_coll.hpp:590
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:3104
int RT_dof[Geometry::NumGeom]
Definition fe_coll.hpp:591
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:3153
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:599
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:602
const char * Name() const override
Definition fe_coll.hpp:608
Arbitrary order 3D H(div)-conforming Raviart-Thomas finite elements in 2D.
Definition fe_coll.hpp:663
void InitFaces(const int p, const int dim, const int map_type, const bool signs)
Definition fe_coll.cpp:3385
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:3371
FiniteElement * RT_Elements[Geometry::NumGeom]
Definition fe_coll.hpp:667
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:686
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.hpp:689
int RT_dof[Geometry::NumGeom]
Definition fe_coll.hpp:668
const char * Name() const override
Definition fe_coll.hpp:695
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:3424
FiniteElementCollection * GetTraceCollection() const override
Definition fe_coll.cpp:3434
int GetContType() const override
Definition fe_coll.hpp:697
Arbitrary order 3D "H^{-1/2}-conforming" face finite elements defined on the interface between mesh e...
Definition fe_coll.hpp:708
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:3460
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
Definition fe_coll.hpp:466
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:2725
FiniteElementCollection * Clone(int p) const override
Instantiate a new collection of the same type with a different order.
Definition fe_coll.hpp:472
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:1381
const FiniteElement * FiniteElementForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1515
int DofForGeometry(Geometry::Type GeomType) const override
Definition fe_coll.cpp:1533
int GetContType() const override
Definition fe_coll.hpp:1403
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:1549
const char * Name() const override
Definition fe_coll.hpp:1401
A 3D refined tri-linear element on a cube.
Class for finite elements with basis functions that return scalar values.
Definition fe_base.hpp:665
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:403
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)