Dune-Functions 2.11
Loading...
Searching...
No Matches
nedelecbasis.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
5// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later
6
7#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
8#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
9
10#include <array>
11#include <dune/common/exceptions.hh>
12
13#include <dune/grid/common/capabilities.hh>
14#include <dune/grid/common/mcmgmapper.hh>
15
16#include <dune/localfunctions/common/localfiniteelementvariant.hh>
17#include <dune/localfunctions/nedelec.hh>
18
23
24namespace Dune::Functions
25{
26
27namespace Impl
28{
29 template<typename GV, int dim, typename R, std::size_t order>
30 class Nedelec1stKindLocalFiniteElementMap
31 {
32 using D = typename GV::ctype;
33 constexpr static bool hasFixedElementType = Capabilities::hasSingleGeometryType<typename GV::Grid>::v;
34
35 using CubeFiniteElement = Nedelec1stKindCubeLocalFiniteElement<D,R,dim,order>;
36 using SimplexFiniteElement = Nedelec1stKindSimplexLocalFiniteElement<D,R,dim,order>;
37
38 public:
39
40 using T = LocalBasisTraits<D, dim, FieldVector<D,dim>, R, dim, FieldVector<R,dim>, FieldMatrix<D,dim,dim> >;
41
42 constexpr static unsigned int topologyId = Capabilities::hasSingleGeometryType<typename GV::Grid>::topologyId; // meaningless if hasFixedElementType is false
43 constexpr static GeometryType type = GeometryType(topologyId, GV::dimension);
44
45 using FiniteElement = std::conditional_t<hasFixedElementType,
46 std::conditional_t<type.isCube(),CubeFiniteElement,SimplexFiniteElement>,
47 LocalFiniteElementVariant<CubeFiniteElement, SimplexFiniteElement> >;
48
49 static std::size_t numVariants(GeometryType type)
50 {
51 if (order!=1) // I am not sure whether the formula below is correct for all orders.
52 DUNE_THROW(NotImplemented, "Only Nedelec elements of order 1 are implemented!");
53
54 auto numEdges = referenceElement<D,dim>(type).size(dim-1);
55 return power(2,numEdges);
56 }
57
58 Nedelec1stKindLocalFiniteElementMap(const GV& gv)
59 : elementMapper_(gv, mcmgElementLayout()),
60 orientation_(gv.size(0))
61 {
62 update(gv);
63 }
64
65 void update(const GV& gv)
66 {
67 elementMapper_.update(gv);
68 orientation_.resize(gv.size(0));
69
70 // create all variants
71 if constexpr (hasFixedElementType)
72 {
73 variants_.resize(numVariants(type));
74 for (size_t i = 0; i < numVariants(type); i++)
75 variants_[i] = FiniteElement(i);
76 }
77 else
78 {
79 // for mixed grids add offset for cubes
80 variants_.resize(numVariants(GeometryTypes::simplex(dim)) + numVariants(GeometryTypes::cube(dim)));
81 for (size_t i = 0; i < numVariants(GeometryTypes::simplex(dim)); i++)
82 variants_[i] = SimplexFiniteElement(i);
83 for (size_t i = 0; i < numVariants(GeometryTypes::cube(dim)); i++)
84 variants_[i + numVariants(GeometryTypes::simplex(dim))] = CubeFiniteElement(i);
85 }
86
87
88 // compute orientation for all elements
89 const auto& indexSet = gv.indexSet();
90
91 for(const auto& element : elements(gv))
92 {
93 const auto& refElement = referenceElement(element);
94 auto elementIndex = elementMapper_.index(element);
95 orientation_[elementIndex] = 0;
96
97 for (std::size_t i=0; i<element.subEntities(dim-1); i++)
98 {
99 // Local vertex indices within the element
100 auto localV0 = refElement.subEntity(i,dim-1, 0,dim);
101 auto localV1 = refElement.subEntity(i,dim-1, 1,dim);
102
103 // Global vertex indices within the grid
104 auto globalV0 = indexSet.subIndex(element,localV0,dim);
105 auto globalV1 = indexSet.subIndex(element,localV1,dim);
106
107 if ( (localV0<localV1 && globalV0>globalV1) || (localV0>localV1 && globalV0<globalV1) )
108 orientation_[elementIndex] |= (1 << i);
109 }
110 // for mixed grids add offset for cubes
111 if constexpr (!hasFixedElementType)
112 if (element.type().isCube())
113 orientation_[elementIndex] += numVariants(GeometryTypes::simplex(dim));
114 }
115 }
116
117 template<class Element>
118 const auto& find(const Element& element) const
119 {
120 return variants_[orientation_[elementMapper_.index(element)]];
121 }
122
123 private:
124 std::vector<FiniteElement> variants_;
125 Dune::MultipleCodimMultipleGeomTypeMapper<GV> elementMapper_;
126 std::vector<unsigned short> orientation_;
127 };
128
129
130} // namespace Impl
131
132
133// *****************************************************************************
134// This is the reusable part of the basis. It contains
135//
136// NedelecPreBasis
137// NedelecNode
138//
139// The pre-basis allows to create the others and is the owner of possible shared
140// state. These components do _not_ depend on the global basis and local view
141// and can be used without a global basis.
142// *****************************************************************************
143
144template<typename GV, typename Range, std::size_t kind, int order>
145class NedelecNode;
146
147template<typename GV, typename Range, std::size_t kind, int order>
149 public LeafPreBasisMapperMixin<GV>
150{
151 static const int dim = GV::dimension;
152 static_assert(kind==1, "Only the Nedelec basis of the first kind is currently implemented!");
153 using FiniteElementMap = typename Impl::Nedelec1stKindLocalFiniteElementMap<GV, dim, Range, order>;
154
155 using Mapper = Dune::MultipleCodimMultipleGeomTypeMapper<GV>;
156public:
157
159 using GridView = GV;
160 using size_type = std::size_t;
161
163
166 : LeafPreBasisMapperMixin<GV>(gv, mcmgLayout(Dim<1>{}))
168 {
169 if (kind!=1)
170 DUNE_THROW(NotImplemented, "Only Nedelec elements of the first kind are implemented!");
171
172 // There is no inherent reason why the basis shouldn't work for grids with more than two
173 // element types. Somebody simply has to sit down and implement the missing bits.
174 if (gv.indexSet().types(0).size() > 2)
175 DUNE_THROW(NotImplemented, "Nédélec basis is only implemented for grids with simplex and cube elements");
176
177 for(auto type : gv.indexSet().types(0))
178 if (!type.isSimplex() && !type.isCube())
179 DUNE_THROW(NotImplemented, "Nédélec basis is only implemented for grids with simplex or cube elements.");
180
181 if (order>1)
182 DUNE_THROW(NotImplemented, "Only first-order elements are implemented");
183
184 if (dim!=2 && dim!=3)
185 DUNE_THROW(NotImplemented, "Only 2d and 3d Nédélec elements are implemented");
186 }
187
189 void update (const GridView& gv)
190 {
192 finiteElementMap_.update(gv);
193 }
194
199 {
200 return Node{&finiteElementMap_};
201 }
202
203protected:
204 FiniteElementMap finiteElementMap_;
205};
206
207
208
209template<typename GV, typename Range, size_t kind, int order>
211 public LeafBasisNode
212{
213 static const int dim = GV::dimension;
214
215public:
216
217 using size_type = std::size_t;
218 using Element = typename GV::template Codim<0>::Entity;
219 static_assert(kind==1, "Only Nedelec elements of the first kind are implemented!");
220 using FiniteElementMap = typename Impl::Nedelec1stKindLocalFiniteElementMap<GV, dim, Range, order>;
221 using FiniteElement = Impl::GlobalValuedLocalFiniteElement<Impl::CovariantPiolaTransformator,
222 typename FiniteElementMap::FiniteElement,
223 Element>;
224
225 NedelecNode(const FiniteElementMap* finiteElementMap) :
226 element_(nullptr),
227 finiteElementMap_(finiteElementMap)
228 { }
229
231 const Element& element() const
232 {
233 return *element_;
234 }
235
241 {
242 return finiteElement_;
243 }
244
246 void bind(const Element& e)
247 {
248 element_ = &e;
249 finiteElement_.bind((finiteElementMap_->find(*element_)), e);
250 this->setSize(finiteElement_.size());
251 }
252
253protected:
254
258};
259
260
261
262namespace BasisFactory {
263
273template<std::size_t kind, std::size_t order, typename Range=double>
275{
276 return [](const auto& gridView) {
277 return NedelecPreBasis<std::decay_t<decltype(gridView)>, Range, kind, order>(gridView);
278 };
279}
280
281} // end namespace BasisFactory
282
283
284
285// *****************************************************************************
286// This is the actual global basis implementation based on the reusable parts.
287// *****************************************************************************
288
296template<typename GV, std::size_t kind, std::size_t order, typename Range=double>
298
299} // end namespace Dune::Functions
300
301
302#endif // DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
auto power(ChildPreBasisFactory &&childPreBasisFactory, std::size_t k, const IndexMergingStrategy &)
Create a pre-basis factory that can build a PowerPreBasis.
Definition dynamicpowerbasis.hh:410
auto nedelec()
Create a pre-basis factory that can create a Nédélec pre-basis.
Definition nedelecbasis.hh:274
auto elements(const SubDomainGridView< HostGridView > &subDomainGridView)
ADL findable access to element range for a SubDomainGridView.
Definition subdomain.hh:487
Definition monomialset.hh:19
DefaultGlobalBasis< NedelecPreBasis< GV, Range, kind, order > > NedelecBasis
Basis of a k-th-order Nédélec finite element space.
Definition nedelecbasis.hh:297
Definition argyrisbasis.hh:926
Global basis for given pre-basis.
Definition defaultglobalbasis.hh:53
void update(const GridView &gv)
Update the stored GridView.
Definition leafprebasismappermixin.hh:101
LeafPreBasisMapperMixin(const GridView &gv, Dune::MCMGLayout layout)
Construct from GridView and local DOF layout.
Definition leafprebasismappermixin.hh:74
Definition nedelecbasis.hh:212
const FiniteElementMap * finiteElementMap_
Definition nedelecbasis.hh:257
FiniteElement finiteElement_
Definition nedelecbasis.hh:255
Impl::GlobalValuedLocalFiniteElement< Impl::CovariantPiolaTransformator, typename FiniteElementMap::FiniteElement, Element > FiniteElement
Definition nedelecbasis.hh:221
void bind(const Element &e)
Bind to element.
Definition nedelecbasis.hh:246
const Element & element() const
Return current element, throw if unbound.
Definition nedelecbasis.hh:231
typename GV::template Codim< 0 >::Entity Element
Definition nedelecbasis.hh:218
const FiniteElement & finiteElement() const
Return the LocalFiniteElement for the element we are bound to.
Definition nedelecbasis.hh:240
NedelecNode(const FiniteElementMap *finiteElementMap)
Definition nedelecbasis.hh:225
typename Impl::Nedelec1stKindLocalFiniteElementMap< GV, dim, Range, order > FiniteElementMap
Definition nedelecbasis.hh:220
const Element * element_
Definition nedelecbasis.hh:256
std::size_t size_type
Definition nedelecbasis.hh:217
Definition nedelecbasis.hh:150
std::size_t size_type
Definition nedelecbasis.hh:160
NedelecPreBasis(const GridView &gv)
Constructor for a given grid view object.
Definition nedelecbasis.hh:165
GV GridView
The grid view that the FE space is defined on.
Definition nedelecbasis.hh:159
FiniteElementMap finiteElementMap_
Definition nedelecbasis.hh:204
NedelecNode< GV, Range, kind, order > Node
Definition nedelecbasis.hh:162
void update(const GridView &gv)
Update the stored grid view, to be called if the grid has changed.
Definition nedelecbasis.hh:189
Node makeNode() const
Create tree node.
Definition nedelecbasis.hh:198
void setSize(const size_type size)
Definition nodes.hh:197
Definition nodes.hh:218