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