DUNE PDELab (git)

qkdg.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_PDELAB_FINITEELEMENTMAP_QKDG_HH
5#define DUNE_PDELAB_FINITEELEMENTMAP_QKDG_HH
6
7#include <dune/pdelab/finiteelement/qkdglagrange.hh>
8#include <dune/pdelab/finiteelement/qkdglegendre.hh>
9#include <dune/pdelab/finiteelement/qkdglobatto.hh>
10
11#include <dune/pdelab/finiteelementmap/opbfem.hh>
12#include <dune/pdelab/finiteelementmap/finiteelementmap.hh>
13
14namespace Dune {
15 namespace PDELab {
16
18 enum class QkDGBasisPolynomial { lagrange, legendre, lobatto, l2orthonormal};
19
20#ifndef DOXYGEN
21 // Class declaration. Use template specialization below.
22 template<class D, class R, int k, int d, QkDGBasisPolynomial p = QkDGBasisPolynomial::lagrange>
23 class QkDGLocalFiniteElementMap;
24#endif
25
36 template<class D, class R, int k, int d>
37 class QkDGLocalFiniteElementMap<D,R,k,d,QkDGBasisPolynomial::lagrange>
38 : public Dune::PDELab::SimpleLocalFiniteElementMap< Dune::QkDGLagrangeLocalFiniteElement<D,R,k,d>,d>
39 {
40 public:
41
42 static constexpr bool fixedSize()
43 {
44 return true;
45 }
46
47 static constexpr bool hasDOFs(int codim)
48 {
49 return codim == 0;
50 }
51
52 static constexpr std::size_t size(GeometryType gt)
53 {
54 if (gt == GeometryTypes::cube(d))
55 return Dune::QkStuff::QkSize<k,d>::value;
56 else
57 return 0;
58 }
59
60 static constexpr std::size_t maxLocalSize()
61 {
62 return Dune::QkStuff::QkSize<k,d>::value;
63 }
64
66 static constexpr QkDGBasisPolynomial polynomial()
67 {
68 return QkDGBasisPolynomial::lagrange;
69 }
70
72 static constexpr std::size_t order()
73 {
74 return k;
75 }
76 };
77
88 template<class D, class R, int k, int d>
89 class QkDGLocalFiniteElementMap<D,R,k,d,QkDGBasisPolynomial::legendre>
90 : public Dune::PDELab::SimpleLocalFiniteElementMap< Dune::QkDGLegendreLocalFiniteElement<D,R,k,d>,d>
91 {
92 public:
93
94 static constexpr bool fixedSize()
95 {
96 return true;
97 }
98
99 static constexpr bool hasDOFs(int codim)
100 {
101 return codim == 0;
102 }
103
104 static constexpr std::size_t size(GeometryType gt)
105 {
106 if (gt == GeometryTypes::cube(d))
107 return Dune::QkStuff::QkSize<k,d>::value;
108 else
109 return 0;
110 }
111
112 static constexpr std::size_t maxLocalSize()
113 {
114 return Dune::LegendreStuff::LegendreSize<k,d>::value;
115 }
116
118 static constexpr QkDGBasisPolynomial polynomial()
119 {
120 return QkDGBasisPolynomial::legendre;
121 }
122
124 static constexpr std::size_t order()
125 {
126 return k;
127 }
128 };
129
140 template<class D, class R, int k, int d>
141 class QkDGLocalFiniteElementMap<D,R,k,d,QkDGBasisPolynomial::lobatto>
142 : public Dune::PDELab::SimpleLocalFiniteElementMap< Dune::QkDGGLLocalFiniteElement<D,R,k,d>,d>
143 {
144 public:
145
146 static constexpr bool fixedSize()
147 {
148 return true;
149 }
150
151 static constexpr bool hasDOFs(int codim)
152 {
153 return codim == 0;
154 }
155
156 static constexpr std::size_t size(GeometryType gt)
157 {
158 if (gt == GeometryTypes::cube(d))
159 return Dune::QkStuff::QkSize<k,d>::value;
160 else
161 return 0;
162 }
163
164 static constexpr std::size_t maxLocalSize()
165 {
166 return Dune::QkStuff::QkSize<k,d>::value;
167 }
168
170 static constexpr QkDGBasisPolynomial polynomial()
171 {
172 return QkDGBasisPolynomial::lobatto;
173 }
174
176 static constexpr std::size_t order()
177 {
178 return k;
179 }
180 };
181
182
199 template<class D, class R, int k, int d>
200 class QkDGLocalFiniteElementMap<D,R,k,d,QkDGBasisPolynomial::l2orthonormal>
201 : public OPBLocalFiniteElementMap<D,R,k,d,Dune::GeometryType::cube,
202#if HAVE_GMP
203 Dune::GMPField<512>,
204#else
205 R,
206#endif
207 Dune::PB::BasisType::Qk>
208 {
209 public:
210
212 static constexpr QkDGBasisPolynomial polynomial()
213 {
214 return QkDGBasisPolynomial::l2orthonormal;
215 }
216
218 static constexpr std::size_t order()
219 {
220 return k;
221 }
222
223 };
224
225 }
226}
227
228#endif // DUNE_PDELAB_FINITEELEMENTMAP_QKDG_HH
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
static constexpr std::size_t order()
return order of polynomial basis
Definition: qkdg.hh:218
static constexpr QkDGBasisPolynomial polynomial()
return type of polynomial basis
Definition: qkdg.hh:212
static constexpr std::size_t order()
return order of polynomial basis
Definition: qkdg.hh:72
static constexpr QkDGBasisPolynomial polynomial()
return type of polynomial basis
Definition: qkdg.hh:66
static constexpr QkDGBasisPolynomial polynomial()
return type of polynomial basis
Definition: qkdg.hh:118
static constexpr std::size_t order()
return order of polynomial basis
Definition: qkdg.hh:124
static constexpr std::size_t order()
return order of polynomial basis
Definition: qkdg.hh:176
static constexpr QkDGBasisPolynomial polynomial()
return type of polynomial basis
Definition: qkdg.hh:170
simple implementation where all entities have the same finite element
Definition: finiteelementmap.hh:101
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
auto lagrange()
Create a pre-basis factory that can create a Lagrange pre-basis.
Definition: lagrangebasis.hh:477
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:462
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
std::size_t fixedSize
The number of data items per index if it is fixed, 0 otherwise.
Definition: variablesizecommunicator.hh:264
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 23, 23:29, 2024)