3#ifndef DUNE_LOCALFUNCTIONS_NEDELEC_NEDELEC1STKINDCUBE_HH
4#define DUNE_LOCALFUNCTIONS_NEDELEC_NEDELEC1STKINDCUBE_HH
12#include <dune/geometry/referenceelements.hh>
15#include <dune/localfunctions/common/localbasis.hh>
16#include <dune/localfunctions/common/localfiniteelementtraits.hh>
17#include <dune/localfunctions/common/localinterpolation.hh>
18#include <dune/localfunctions/common/localkey.hh>
34 template<
class D,
class R,
int dim,
int k>
35 class Nedelec1stKindCubeLocalBasis
38 constexpr static std::size_t numberOfEdges =
power(2,dim-1)*dim;
41 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,
42 R,dim,FieldVector<R,dim>,
43 FieldMatrix<R,dim,dim> >;
51 Nedelec1stKindCubeLocalBasis()
53 std::fill(edgeOrientation_.begin(), edgeOrientation_.end(), 1.0);
58 Nedelec1stKindCubeLocalBasis(std::bitset<numberOfEdges> edgeOrientation)
59 : Nedelec1stKindCubeLocalBasis()
61 for (std::size_t i=0; i<edgeOrientation_.size(); i++)
62 edgeOrientation_[i] *= edgeOrientation[i] ? -1.0 : 1.0;
66 static constexpr unsigned int size()
68 static_assert(dim==2 || dim==3,
"Nedelec shape functions are implemented only for 2d and 3d cubes.");
72 return 3*k * (k+1) * (k+1);
81 std::vector<typename Traits::RangeType>& out)
const
83 static_assert(k==1,
"Evaluating Nédélec shape functions is implemented only for first order.");
100 out[0] = { 0, D(1) - in[0]};
101 out[1] = { 0, in[0]};
102 out[2] = { D(1) - in[1], 0};
103 out[3] = { in[1], 0};
106 if constexpr (dim==3)
125 for (std::size_t i=0; i<out.size(); i++)
128 out[0][2] = { 1 - in[0] - in[1] + in[0]*in[1]};
129 out[1][2] = { in[0] - in[0]*in[1]};
130 out[2][2] = { in[1] - in[0]*in[1]};
131 out[3][2] = { in[0]*in[1]};
133 out[4][1] = { 1 - in[0] - in[2] + in[0]*in[2]};
134 out[5][1] = { in[0] - in[0]*in[2]};
135 out[8][1] = { in[2] - in[0]*in[2]};
136 out[9][1] = { in[0]*in[2]};
138 out[6][0] = { 1 - in[1] - in[2] + in[1]*in[2]};
139 out[7][0] = { in[1] - in[1]*in[2]};
140 out[10][0] = { in[2] - in[1]*in[2]};
141 out[11][0] = { in[1]*in[2]};
144 for (std::size_t i=0; i<out.size(); i++)
145 out[i] *= edgeOrientation_[i];
154 std::vector<typename Traits::JacobianType>& out)
const
159 for (std::size_t i=0; i<out.size(); i++)
160 for (std::size_t j=0; j<dim; j++)
163 out[0][1] = { -1, 0};
166 out[2][0] = { 0, -1};
171 for (std::size_t i=0; i<out.size(); i++)
172 for(std::size_t j=0;j<dim; j++)
176 out[0][2] = {-1 +in[1], -1 + in[0], 0};
177 out[1][2] = { 1 -in[1], - in[0], 0};
178 out[2][2] = { -in[1], 1 - in[0], 0};
179 out[3][2] = { in[1], in[0], 0};
181 out[4][1] = {-1 +in[2], 0, -1 + in[0]};
182 out[5][1] = { 1 -in[2], 0, - in[0]};
183 out[8][1] = { -in[2], 0, 1 - in[0]};
184 out[9][1] = { in[2], 0, in[0]};
186 out[6][0] = { 0, -1 + in[2], -1 + in[1]};
187 out[7][0] = { 0, 1 - in[2], - in[1]};
188 out[10][0] = { 0, - in[2], 1 - in[1]};
189 out[11][0] = { 0, in[2], in[1]};
193 for (std::size_t i=0; i<out.size(); i++)
194 out[i] *= edgeOrientation_[i];
204 void partial(
const std::array<unsigned int, dim>& order,
206 std::vector<typename Traits::RangeType>& out)
const
209 if (totalOrder == 0) {
210 evaluateFunction(in, out);
211 }
else if (totalOrder == 1) {
212 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
238 out[0] = { 0, 0, -1 +in[1]};
239 out[1] = { 0, 0, 1 -in[1]};
240 out[2] = { 0, 0, -in[1]};
241 out[3] = { 0, 0, in[1]};
243 out[4] = { 0, -1 +in[2], 0};
244 out[5] = { 0, 1 -in[2], 0};
245 out[8] = { 0, -in[2], 0};
246 out[9] = { 0, in[2], 0};
255 out[0] = { 0, 0, -1 + in[0]};
256 out[1] = { 0, 0, - in[0]};
257 out[2] = { 0, 0, 1 - in[0]};
258 out[3] = { 0, 0, in[0]};
265 out[6] = { -1 + in[2], 0, 0};
266 out[7] = { 1 - in[2], 0, 0};
267 out[10] = { - in[2], 0, 0};
268 out[11] = { in[2], 0, 0};
277 out[4] = { 0, -1 + in[0], 0};
278 out[5] = { 0, - in[0], 0};
279 out[8] = { 0, 1 - in[0], 0};
280 out[9] = { 0, in[0], 0};
282 out[6] = { -1 + in[1], 0, 0};
283 out[7] = { - in[1], 0, 0};
284 out[10] = { 1 - in[1], 0, 0};
285 out[11] = { in[1], 0, 0};
290 for (std::size_t i=0; i<out.size(); i++)
291 out[i] *= edgeOrientation_[i];
293 }
else if (totalOrder == 2) {
297 for (std::size_t i = 0; i < size(); ++i)
298 for (std::size_t j = 0; j < dim; ++j)
303 for(
size_t i=0; i<out.size(); i++)
307 if( order[0] == 1 and order[1]==1)
310 out[1] = { 0, 0, -1};
311 out[2] = { 0, 0, -1};
316 if( order[0] == 1 and order[2]==1)
319 out[5] = { 0, -1, 0};
320 out[8] = { 0, -1, 0};
325 if( order[1] == 1 and order[2]==1)
328 out[7] = { -1, 0, 0};
329 out[10] = { -1, 0, 0};
330 out[11] = { 1, 0, 0};
333 for (std::size_t i=0; i<out.size(); i++)
334 out[i] *= edgeOrientation_[i];
340 for (std::size_t i = 0; i < size(); ++i)
341 for (std::size_t j = 0; j < dim; ++j)
348 unsigned int order()
const
359 std::array<R,numberOfEdges> edgeOrientation_;
367 template <
int dim,
int k>
368 class Nedelec1stKindCubeLocalCoefficients
372 Nedelec1stKindCubeLocalCoefficients ()
375 static_assert(k==1,
"Only first-order Nédélec local coefficients are implemented.");
378 for (std::size_t i=0; i<size(); i++)
379 localKey_[i] = LocalKey(i,dim-1,0);
383 std::size_t size()
const
385 static_assert(dim==2 || dim==3,
"Nédélec shape functions are implemented only for 2d and 3d cubes.");
386 return (dim==2) ? 2*k * (k+1)
387 : 3*k * (k+1) * (k+1);
392 const LocalKey& localKey (std::size_t i)
const
398 std::vector<LocalKey> localKey_;
406 class Nedelec1stKindCubeLocalInterpolation
408 static constexpr auto dim = LB::Traits::dimDomain;
409 static constexpr auto size = LB::size();
412 constexpr static std::size_t numberOfEdges =
power(2,dim-1)*dim;
417 Nedelec1stKindCubeLocalInterpolation (std::bitset<numberOfEdges> s = 0)
421 for (std::size_t i=0; i<numberOfEdges; i++)
422 m_[i] = refElement.position(i,dim-1);
424 for (std::size_t i=0; i<numberOfEdges; i++)
426 auto vertexIterator = refElement.subEntities(i,dim-1,dim).begin();
427 auto v0 = *vertexIterator;
428 auto v1 = *(++vertexIterator);
434 edge_[i] = refElement.position(v1,dim) - refElement.position(v0,dim);
435 edge_[i] *= (s[i]) ? -1.0 : 1.0;
444 template<
typename F,
typename C>
445 void interpolate (
const F& ff, std::vector<C>& out)
const
448 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
450 for (std::size_t i=0; i<size; i++)
454 for (
int j=0; j<dim; j++)
455 out[i] += y[j]*edge_[i][j];
461 std::array<typename LB::Traits::DomainType, numberOfEdges> m_;
463 std::array<typename LB::Traits::DomainType, numberOfEdges> edge_;
492 template<
class D,
class R,
int dim,
int k>
497 Impl::Nedelec1stKindCubeLocalCoefficients<dim,k>,
498 Impl::Nedelec1stKindCubeLocalInterpolation<Impl::Nedelec1stKindCubeLocalBasis<D,R,dim,k> > >;
500 static_assert(dim==2 || dim==3,
"Nedelec elements are only implemented for 2d and 3d elements.");
501 static_assert(k==1,
"Nedelec elements of the first kind are currently only implemented for order k==1.");
524 return coefficients_;
529 return interpolation_;
532 static constexpr unsigned int size ()
534 return Traits::LocalBasisType::size();
Nédélec elements of the first kind for cube elements.
Definition: nedelec1stkindcube.hh:494
Nedelec1stKindCubeLocalFiniteElement()=default
Default constructor.
Nedelec1stKindCubeLocalFiniteElement(std::bitset< power(2, dim-1) *dim > s)
Constructor with explicitly given edge orientations.
Definition: nedelec1stkindcube.hh:512
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:130
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:470
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:289
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:11
constexpr Mantissa power(Mantissa m, Exponent p)
Power method for integer exponents.
Definition: math.hh:73
D DomainType
domain type
Definition: localbasis.hh:43
traits helper struct
Definition: localfiniteelementtraits.hh:11
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
A unique label for each type of element that can occur in a grid.