3#ifndef DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
4#define DUNE_FUNCTIONS_FUNCTIONSPACEBASES_NEDELECBASIS_HH
12#include <dune/localfunctions/common/localfiniteelementvariant.hh>
13#include <dune/localfunctions/nedelec.hh>
15#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
16#include <dune/functions/functionspacebases/flatmultiindex.hh>
17#include <dune/functions/functionspacebases/globalvaluedlocalfiniteelement.hh>
18#include <dune/functions/functionspacebases/nodes.hh>
20namespace Dune::Functions
25 template<
typename GV,
int dim,
typename R, std::
size_t order>
26 class Nedelec1stKindLocalFiniteElementMap
28 using D =
typename GV::ctype;
29 constexpr static bool hasFixedElementType = Capabilities::hasSingleGeometryType<typename GV::Grid>::v;
31 using CubeFiniteElement = Nedelec1stKindCubeLocalFiniteElement<D,R,dim,order>;
32 using SimplexFiniteElement = Nedelec1stKindSimplexLocalFiniteElement<D,R,dim,order>;
36 using T = LocalBasisTraits<D, dim, FieldVector<D,dim>, R, dim, FieldVector<R,dim>, FieldMatrix<D,dim,dim> >;
38 constexpr static unsigned int topologyId = Capabilities::hasSingleGeometryType<typename GV::Grid>::topologyId;
41 using FiniteElement = std::conditional_t<hasFixedElementType,
42 std::conditional_t<type.isCube(),CubeFiniteElement,SimplexFiniteElement>,
43 LocalFiniteElementVariant<CubeFiniteElement, SimplexFiniteElement> >;
45 static std::size_t numVariants(GeometryType type)
48 DUNE_THROW(NotImplemented,
"Only Nedelec elements of order 1 are implemented!");
50 auto numEdges = referenceElement<D,dim>(type).size(dim-1);
51 return power(2,numEdges);
54 Nedelec1stKindLocalFiniteElementMap(
const GV& gv)
56 orientation_(gv.size(0))
59 if constexpr (hasFixedElementType)
61 variants_.resize(numVariants(type));
62 for (
size_t i = 0; i < numVariants(type); i++)
63 variants_[i] = FiniteElement(i);
70 variants_[i] = SimplexFiniteElement(i);
77 const auto& indexSet = gv.indexSet();
79 for(
const auto& element : elements(gv))
82 auto elementIndex = elementMapper_.index(element);
83 orientation_[elementIndex] = 0;
85 for (std::size_t i=0; i<element.subEntities(dim-1); i++)
88 auto localV0 = refElement.subEntity(i,dim-1, 0,dim);
89 auto localV1 = refElement.subEntity(i,dim-1, 1,dim);
92 auto globalV0 = indexSet.subIndex(element,localV0,dim);
93 auto globalV1 = indexSet.subIndex(element,localV1,dim);
95 if ( (localV0<localV1 && globalV0>globalV1) || (localV0>localV1 && globalV0<globalV1) )
96 orientation_[elementIndex] |= (1 << i);
99 if constexpr (!hasFixedElementType)
100 if (element.type().isCube())
105 template<
class Element>
106 const auto& find(
const Element& element)
const
108 return variants_[orientation_[elementMapper_.index(element)]];
112 std::vector<FiniteElement> variants_;
114 std::vector<unsigned short> orientation_;
132template<
typename GV,
typename Range, std::
size_t kind,
int order>
135template<
typename GV,
typename Range, std::
size_t kind,
int order,
class MI>
138 static const int dim = GV::dimension;
139 static_assert(kind==1,
"Only the Nedelec basis of the first kind is currently implemented!");
140 using FiniteElementMap =
typename Impl::Nedelec1stKindLocalFiniteElementMap<GV, dim, Range, order>;
147 using size_type = std::size_t;
149 using Node = NedelecNode<GV, Range, kind, order>;
152 using IndexSet = Impl::DefaultNodeIndexSet<NedelecPreBasis>;
155 using MultiIndex = MI;
160 NedelecPreBasis(
const GridView& gv) :
162 finiteElementMap_(gv),
166 DUNE_THROW(NotImplemented,
"Only Nedelec elements of the first kind are implemented!");
170 if (gv.indexSet().types(0).size() > 2)
171 DUNE_THROW(NotImplemented,
"Nédélec basis is only implemented for grids with simplex and cube elements");
173 for(
auto type : gv.indexSet().types(0))
174 if (!type.isSimplex() && !type.isCube())
175 DUNE_THROW(NotImplemented,
"Nédélec basis is only implemented for grids with simplex or cube elements.");
178 DUNE_THROW(NotImplemented,
"Only first-order elements are implemented");
180 if (dim!=2 && dim!=3)
181 DUNE_THROW(NotImplemented,
"Only 2d and 3d Nédélec elements are implemented");
184 void initializeIndices()
189 const GridView& gridView()
const
195 void update (
const GridView& gv)
198 mapper_.update(gridView_);
204 Node makeNode()
const
206 return Node{&finiteElementMap_};
216 [[deprecated(
"Warning: The IndexSet typedef and the makeIndexSet method are deprecated. "\
217 "As a replacement use the indices() method of the PreBasis directly.")]]
218 IndexSet makeIndexSet()
const
220 return IndexSet{*
this};
223 size_type size()
const
225 return mapper_.size();
229 size_type size(
const SizePrefix prefix)
const
231 assert(prefix.size() == 0 || prefix.size() == 1);
232 return (prefix.size() == 0) ? size() : 0;
235 size_type dimension()
const
240 size_type maxNodeSize()
const
242 size_type result = 0;
243 for (
auto&& type : gridView_.indexSet().types(0))
245 size_type numEdges = referenceElement<typename GV::ctype,dim>(type).size(dim-1);
246 result =
std::max(result, numEdges);
255 template<
typename It>
256 It indices(
const Node& node, It it)
const
258 const auto& element = node.element();
261 if (not(element.type().isCube()) and not(element.type().isSimplex()))
262 DUNE_THROW(NotImplemented,
"NedelecBasis only implemented for cube and simplex elements.");
264 for(std::size_t i=0, end=node.size(); i<end; ++i, ++it)
266 Dune::LocalKey localKey = node.finiteElement().localCoefficients().localKey(i);
267 *it = { mapper_.subIndex(element, localKey.
subEntity(), localKey.
codim()) + localKey.
index() };
275 FiniteElementMap finiteElementMap_;
281template<
typename GV,
typename Range,
size_t kind,
int order>
285 static const int dim = GV::dimension;
289 using size_type = std::size_t;
290 using Element =
typename GV::template Codim<0>::Entity;
291 static_assert(kind==1,
"Only Nedelec elements of the first kind are implemented!");
292 using FiniteElementMap =
typename Impl::Nedelec1stKindLocalFiniteElementMap<GV, dim, Range, order>;
293 using FiniteElement = Impl::GlobalValuedLocalFiniteElement<Impl::CovariantPiolaTransformator,
294 typename FiniteElementMap::FiniteElement,
297 NedelecNode(
const FiniteElementMap* finiteElementMap) :
299 finiteElementMap_(finiteElementMap)
303 const Element& element()
const
312 const FiniteElement& finiteElement()
const
314 return finiteElement_;
318 void bind(
const Element& e)
321 finiteElement_.bind((finiteElementMap_->find(*element_)), e);
322 this->setSize(finiteElement_.size());
327 FiniteElement finiteElement_;
328 const Element* element_;
329 const FiniteElementMap* finiteElementMap_;
334namespace BasisFactory {
338template<std::
size_t kind, std::
size_t order,
typename Range>
339class NedelecPreBasisFactory
342 static const std::size_t requiredMultiIndexSize=1;
344 template<
class MultiIndex,
class Gr
idView>
345 auto makePreBasis(
const GridView& gridView)
const
347 return NedelecPreBasis<GridView, Range, kind, order, MultiIndex>(gridView);
363template<std::
size_t kind, std::
size_t order,
typename Range=
double>
366 return Impl::NedelecPreBasisFactory<kind, order, Range>();
384template<
typename GV, std::
size_t kind, std::
size_t order,
typename Range=
double>
385using NedelecBasis = DefaultGlobalBasis<NedelecPreBasis<GV, Range, kind, order, FlatMultiIndex<std::size_t> > >;
Describe position of one degree of freedom.
Definition: localkey.hh:21
unsigned int index() const
Return offset within subentity.
Definition: localkey.hh:66
unsigned int codim() const
Return codim of associated entity.
Definition: localkey.hh:60
unsigned int subEntity() const
Return number of associated subentity.
Definition: localkey.hh:54
Implementation class for a multiple codim and multiple geometry type mapper.
Definition: mcmgmapper.hh:127
A Vector class with statically reserved memory.
Definition: reservedvector.hh:43
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:130
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
auto nedelec()
Create a pre-basis factory that can create a Nédélec pre-basis.
Definition: nedelecbasis.hh:364
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:470
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:461
MCMGLayout mcmgLayout(Codim< codim >)
layout for entities of codimension codim
Definition: mcmgmapper.hh:70
MCMGLayout mcmgElementLayout()
layout for elements (codim-0 entities)
Definition: mcmgmapper.hh:95
auto max(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::max()
Definition: defaults.hh:79
Mapper for multiple codim and multiple geometry types.
constexpr Mantissa power(Mantissa m, Exponent p)
Power method for integer exponents.
Definition: math.hh:73