DUNE PDELab (git)

nedelecbasis.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
5
6#include <array>
8
11
12#include <dune/localfunctions/common/localfiniteelementvariant.hh>
13#include <dune/localfunctions/nedelec.hh>
14
15#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
16#include <dune/functions/functionspacebases/globalvaluedlocalfiniteelement.hh>
17#include <dune/functions/functionspacebases/leafprebasismixin.hh>
18#include <dune/functions/functionspacebases/nodes.hh>
19
20namespace Dune::Functions
21{
22
23namespace Impl
24{
25 template<typename GV, int dim, typename R, std::size_t order>
26 class Nedelec1stKindLocalFiniteElementMap
27 {
28 using D = typename GV::ctype;
29 constexpr static bool hasFixedElementType = Capabilities::hasSingleGeometryType<typename GV::Grid>::v;
30
31 using CubeFiniteElement = Nedelec1stKindCubeLocalFiniteElement<D,R,dim,order>;
32 using SimplexFiniteElement = Nedelec1stKindSimplexLocalFiniteElement<D,R,dim,order>;
33
34 public:
35
36 using T = LocalBasisTraits<D, dim, FieldVector<D,dim>, R, dim, FieldVector<R,dim>, FieldMatrix<D,dim,dim> >;
37
38 constexpr static unsigned int topologyId = Capabilities::hasSingleGeometryType<typename GV::Grid>::topologyId; // meaningless if hasFixedElementType is false
39 constexpr static GeometryType type = GeometryType(topologyId, GV::dimension);
40
41 using FiniteElement = std::conditional_t<hasFixedElementType,
42 std::conditional_t<type.isCube(),CubeFiniteElement,SimplexFiniteElement>,
43 LocalFiniteElementVariant<CubeFiniteElement, SimplexFiniteElement> >;
44
45 static std::size_t numVariants(GeometryType type)
46 {
47 if (order!=1) // I am not sure whether the formula below is correct for all orders.
48 DUNE_THROW(NotImplemented, "Only Nedelec elements of order 1 are implemented!");
49
50 auto numEdges = referenceElement<D,dim>(type).size(dim-1);
51 return power(2,numEdges);
52 }
53
54 Nedelec1stKindLocalFiniteElementMap(const GV& gv)
55 : elementMapper_(gv, mcmgElementLayout()),
56 orientation_(gv.size(0))
57 {
58 // create all variants
59 if constexpr (hasFixedElementType)
60 {
61 variants_.resize(numVariants(type));
62 for (size_t i = 0; i < numVariants(type); i++)
63 variants_[i] = FiniteElement(i);
64 }
65 else
66 {
67 // for mixed grids add offset for cubes
68 variants_.resize(numVariants(GeometryTypes::simplex(dim)) + numVariants(GeometryTypes::cube(dim)));
69 for (size_t i = 0; i < numVariants(GeometryTypes::simplex(dim)); i++)
70 variants_[i] = SimplexFiniteElement(i);
71 for (size_t i = 0; i < numVariants(GeometryTypes::cube(dim)); i++)
72 variants_[i + numVariants(GeometryTypes::simplex(dim))] = CubeFiniteElement(i);
73 }
74
75
76 // compute orientation for all elements
77 const auto& indexSet = gv.indexSet();
78
79 for(const auto& element : elements(gv))
80 {
81 const auto& refElement = referenceElement(element);
82 auto elementIndex = elementMapper_.index(element);
83 orientation_[elementIndex] = 0;
84
85 for (std::size_t i=0; i<element.subEntities(dim-1); i++)
86 {
87 // Local vertex indices within the element
88 auto localV0 = refElement.subEntity(i,dim-1, 0,dim);
89 auto localV1 = refElement.subEntity(i,dim-1, 1,dim);
90
91 // Global vertex indices within the grid
92 auto globalV0 = indexSet.subIndex(element,localV0,dim);
93 auto globalV1 = indexSet.subIndex(element,localV1,dim);
94
95 if ( (localV0<localV1 && globalV0>globalV1) || (localV0>localV1 && globalV0<globalV1) )
96 orientation_[elementIndex] |= (1 << i);
97 }
98 // for mixed grids add offset for cubes
99 if constexpr (!hasFixedElementType)
100 if (element.type().isCube())
101 orientation_[elementIndex] += numVariants(GeometryTypes::simplex(dim));
102 }
103 }
104
105 template<class Element>
106 const auto& find(const Element& element) const
107 {
108 return variants_[orientation_[elementMapper_.index(element)]];
109 }
110
111 private:
112 std::vector<FiniteElement> variants_;
114 std::vector<unsigned short> orientation_;
115 };
116
117
118} // namespace Impl
119
120
121// *****************************************************************************
122// This is the reusable part of the basis. It contains
123//
124// NedelecPreBasis
125// NedelecNode
126//
127// The pre-basis allows to create the others and is the owner of possible shared
128// state. These components do _not_ depend on the global basis and local view
129// and can be used without a global basis.
130// *****************************************************************************
131
132template<typename GV, typename Range, std::size_t kind, int order>
133class NedelecNode;
134
135template<typename GV, typename Range, std::size_t kind, int order>
136class NedelecPreBasis :
137 public LeafPreBasisMixin< NedelecPreBasis<GV,Range,kind,order> >
138{
139 static const int dim = GV::dimension;
140 static_assert(kind==1, "Only the Nedelec basis of the first kind is currently implemented!");
141 using FiniteElementMap = typename Impl::Nedelec1stKindLocalFiniteElementMap<GV, dim, Range, order>;
142
144public:
145
147 using GridView = GV;
148 using size_type = std::size_t;
149
150 using Node = NedelecNode<GV, Range, kind, order>;
151
153 NedelecPreBasis(const GridView& gv) :
154 gridView_(gv),
155 finiteElementMap_(gv),
156 mapper_(gridView_, mcmgLayout(Dim<1>{}))
157 {
158 if (kind!=1)
159 DUNE_THROW(NotImplemented, "Only Nedelec elements of the first kind are implemented!");
160
161 // There is no inherent reason why the basis shouldn't work for grids with more than two
162 // element types. Somebody simply has to sit down and implement the missing bits.
163 if (gv.indexSet().types(0).size() > 2)
164 DUNE_THROW(NotImplemented, "Nédélec basis is only implemented for grids with simplex and cube elements");
165
166 for(auto type : gv.indexSet().types(0))
167 if (!type.isSimplex() && !type.isCube())
168 DUNE_THROW(NotImplemented, "Nédélec basis is only implemented for grids with simplex or cube elements.");
169
170 if (order>1)
171 DUNE_THROW(NotImplemented, "Only first-order elements are implemented");
172
173 if (dim!=2 && dim!=3)
174 DUNE_THROW(NotImplemented, "Only 2d and 3d Nédélec elements are implemented");
175 }
176
177 void initializeIndices()
178 {}
179
182 const GridView& gridView() const
183 {
184 return gridView_;
185 }
186
187 /* \brief Update the stored grid view, to be called if the grid has changed */
188 void update (const GridView& gv)
189 {
190 gridView_ = gv;
191 mapper_.update(gridView_);
192 }
193
197 Node makeNode() const
198 {
199 return Node{&finiteElementMap_};
200 }
201
202 size_type dimension() const
203 {
204 return mapper_.size();
205 }
206
207 size_type maxNodeSize() const
208 {
209 size_type result = 0;
210 for (auto&& type : gridView_.indexSet().types(0))
211 {
212 size_type numEdges = referenceElement<typename GV::ctype,dim>(type).size(dim-1);
213 result = std::max(result, numEdges);
214 }
215
216 return result;
217 }
218
222 template<typename It>
223 It indices(const Node& node, It it) const
224 {
225 const auto& element = node.element();
226
227 // throw if Element is not of predefined type
228 if (not(element.type().isCube()) and not(element.type().isSimplex()))
229 DUNE_THROW(NotImplemented, "NedelecBasis only implemented for cube and simplex elements.");
230
231 for(std::size_t i=0, end=node.size(); i<end; ++i, ++it)
232 {
233 Dune::LocalKey localKey = node.finiteElement().localCoefficients().localKey(i);
234 *it = { mapper_.subIndex(element, localKey.subEntity(), localKey.codim()) + localKey.index() };
235 }
236
237 return it;
238 }
239
240protected:
241 GridView gridView_;
242 FiniteElementMap finiteElementMap_;
243 Mapper mapper_;
244};
245
246
247
248template<typename GV, typename Range, size_t kind, int order>
249class NedelecNode :
250 public LeafBasisNode
251{
252 static const int dim = GV::dimension;
253
254public:
255
256 using size_type = std::size_t;
257 using Element = typename GV::template Codim<0>::Entity;
258 static_assert(kind==1, "Only Nedelec elements of the first kind are implemented!");
259 using FiniteElementMap = typename Impl::Nedelec1stKindLocalFiniteElementMap<GV, dim, Range, order>;
260 using FiniteElement = Impl::GlobalValuedLocalFiniteElement<Impl::CovariantPiolaTransformator,
261 typename FiniteElementMap::FiniteElement,
262 Element>;
263
264 NedelecNode(const FiniteElementMap* finiteElementMap) :
265 element_(nullptr),
266 finiteElementMap_(finiteElementMap)
267 { }
268
270 const Element& element() const
271 {
272 return *element_;
273 }
274
279 const FiniteElement& finiteElement() const
280 {
281 return finiteElement_;
282 }
283
285 void bind(const Element& e)
286 {
287 element_ = &e;
288 finiteElement_.bind((finiteElementMap_->find(*element_)), e);
289 this->setSize(finiteElement_.size());
290 }
291
292protected:
293
294 FiniteElement finiteElement_;
295 const Element* element_;
296 const FiniteElementMap* finiteElementMap_;
297};
298
299
300
301namespace BasisFactory {
302
312template<std::size_t kind, std::size_t order, typename Range=double>
314{
315 return [](const auto& gridView) {
316 return NedelecPreBasis<std::decay_t<decltype(gridView)>, Range, kind, order>(gridView);
317 };
318}
319
320} // end namespace BasisFactory
321
322
323
324// *****************************************************************************
325// This is the actual global basis implementation based on the reusable parts.
326// *****************************************************************************
327
335template<typename GV, std::size_t kind, std::size_t order, typename Range=double>
336using NedelecBasis = DefaultGlobalBasis<NedelecPreBasis<GV, Range, kind, order > >;
337
338} // end namespace Dune::Functions
339
340
341#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
Describe position of one degree of freedom.
Definition: localkey.hh:24
constexpr unsigned int index() const noexcept
Return offset within subentity.
Definition: localkey.hh:70
constexpr unsigned int codim() const noexcept
Return codim of associated entity.
Definition: localkey.hh:63
constexpr unsigned int subEntity() const noexcept
Return number of associated subentity.
Definition: localkey.hh:56
Implementation class for a multiple codim and multiple geometry type mapper.
Definition: mcmgmapper.hh:129
A set of traits classes to store static information about grid implementation.
A few common exception classes.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
auto nedelec()
Create a pre-basis factory that can create a Nédélec pre-basis.
Definition: nedelecbasis.hh:313
unspecified value type referenceElement(T &&... t)
Returns a reference element for the objects t....
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:462
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:453
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
MCMGLayout mcmgLayout(Codim< codim >)
layout for entities of codimension codim
Definition: mcmgmapper.hh:72
MCMGLayout mcmgElementLayout()
layout for elements (codim-0 entities)
Definition: mcmgmapper.hh:97
Mapper for multiple codim and multiple geometry types.
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
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)