5#ifndef DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
15#include <dune/geometry/referenceelements.hh>
17#include <dune/localfunctions/common/localbasis.hh>
18#include <dune/localfunctions/common/localfiniteelementtraits.hh>
19#include <dune/localfunctions/common/localinterpolation.hh>
20#include <dune/localfunctions/common/localkey.hh>
22namespace Dune {
namespace Impl
25 template<
class LocalBasis>
26 class LagrangeCubeLocalInterpolation;
38 template<
class D,
class R,
unsigned int dim,
unsigned int k>
39 class LagrangeCubeLocalBasis
41 friend class LagrangeCubeLocalInterpolation<LagrangeCubeLocalBasis<D,R,dim,k> >;
44 static R p(
unsigned int i, D x)
47 for (
unsigned int j=0; j<=k; j++)
48 if (j!=i) result *= (k*x-j)/((
int)i-(int)j);
53 static R dp(
unsigned int i, D x)
57 for (
unsigned int j=0; j<=k; j++)
61 R prod( (k*1.0)/((
int)i-(
int)j) );
62 for (
unsigned int l=0; l<=k; l++)
64 prod *= (k*x-l)/((
int)i-(int)l);
73 static R ddp(
unsigned int j, D x)
77 for (
unsigned int i=0; i<=k; i++)
84 for (
unsigned int m=0; m<=k; m++)
89 R prod( (k*1.0)/((
int)j-(
int)m) );
90 for (
unsigned int l=0; l<=k; l++)
91 if (l!=i && l!=j && l!=m)
92 prod *= (k*x-l)/((
int)j-(int)l);
96 result += sum * ( (k*1.0)/((
int)j-(int)i) );
103 static std::array<unsigned int,dim> multiindex (
unsigned int i)
105 std::array<unsigned int,dim> alpha;
106 for (
unsigned int j=0; j<dim; j++)
108 alpha[j] = i % (k+1);
115 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
119 static constexpr unsigned int size ()
121 return power(k+1, dim);
126 std::vector<typename Traits::RangeType>& out)
const
139 for (
size_t i=0; i<size(); i++)
143 for (
unsigned int j=0; j<dim; j++)
145 out[i] *= (i & (1<<j)) ? x[j] : 1-x[j];
151 for (
size_t i=0; i<size(); i++)
154 std::array<unsigned int,dim> alpha(multiindex(i));
160 for (
unsigned int j=0; j<dim; j++)
161 out[i] *= p(alpha[j],x[j]);
171 std::vector<typename Traits::JacobianType>& out)
const
178 std::fill(out[0][0].begin(), out[0][0].end(), 0);
186 for (
size_t i=0; i<size(); i++)
189 for (
unsigned int j=0; j<dim; j++)
193 out[i][0][j] = (i & (1<<j)) ? 1 : -1;
195 for (
unsigned int l=0; l<dim; l++)
199 out[i][0][j] *= (i & (1<<l)) ? x[l] : 1-x[l];
209 for (
size_t i=0; i<size(); i++)
212 std::array<unsigned int,dim> alpha(multiindex(i));
215 for (
unsigned int j=0; j<dim; j++)
219 out[i][0][j] = dp(alpha[j],x[j]);
222 for (
unsigned int l=0; l<dim; l++)
224 out[i][0][j] *= p(alpha[l],x[l]);
235 void partial(
const std::array<unsigned int,dim>& order,
237 std::vector<typename Traits::RangeType>& out)
const
245 out[0] = (totalOrder==0);
253 evaluateFunction(in, out);
255 else if (totalOrder == 1)
259 auto direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
260 if (direction >= dim)
261 DUNE_THROW(RangeError,
"Direction of partial derivative not found!");
264 for (std::size_t i = 0; i < size(); ++i)
268 out[i] = (i & (1<<direction)) ? 1 : -1;
270 for (
unsigned int j = 0; j < dim; ++j)
274 out[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
278 else if (totalOrder == 2)
281 for (
size_t i=0; i<size(); i++)
284 std::array<unsigned int,dim> alpha(multiindex(i));
290 for (std::size_t l=0; l<dim; l++)
295 out[i][0] *= p(alpha[l],in[l]);
299 out[i][0] *= dp(alpha[l],in[l]);
305 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
311 DUNE_THROW(NotImplemented,
"Partial derivative of order " << totalOrder <<
" is not implemented!");
319 for (
size_t i=0; i<size(); i++)
322 std::array<unsigned int,dim> alpha(multiindex(i));
328 for (std::size_t l=0; l<dim; l++)
333 out[i][0] *= p(alpha[l],in[l]);
336 out[i][0] *= dp(alpha[l],in[l]);
339 out[i][0] *= ddp(alpha[l],in[l]);
342 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
349 static constexpr unsigned int order ()
360 template<
unsigned int dim,
unsigned int k>
361 class LagrangeCubeLocalCoefficients
364 static std::array<unsigned int,dim> multiindex (
unsigned int i)
366 std::array<unsigned int,dim> alpha;
367 for (
unsigned int j=0; j<dim; j++)
369 alpha[j] = i % (k+1);
376 void setup1d(std::vector<unsigned int>& subEntity)
380 unsigned lastIndex=0;
387 subEntity[lastIndex++] = 0;
388 for (
unsigned i = 0; i < k - 1; ++i)
389 subEntity[lastIndex++] = 0;
391 subEntity[lastIndex++] = 1;
393 assert(
power(k+1,dim)==lastIndex);
396 void setup2d(std::vector<unsigned int>& subEntity)
400 unsigned lastIndex=0;
414 subEntity[lastIndex++] = 0;
415 for (
unsigned i = 0; i < k - 1; ++i)
416 subEntity[lastIndex++] = 2;
418 subEntity[lastIndex++] = 1;
421 for (
unsigned e = 0; e < k - 1; ++e) {
422 subEntity[lastIndex++] = 0;
423 for (
unsigned i = 0; i < k - 1; ++i)
424 subEntity[lastIndex++] = 0;
425 subEntity[lastIndex++] = 1;
429 subEntity[lastIndex++] = 2;
430 for (
unsigned i = 0; i < k - 1; ++i)
431 subEntity[lastIndex++] = 3;
433 subEntity[lastIndex++] = 3;
435 assert(
power(k+1,dim)==lastIndex);
438 void setup3d(std::vector<unsigned int>& subEntity)
442 unsigned lastIndex=0;
444 const unsigned numIndices =
power(k+1,dim);
445 const unsigned numFaceIndices =
power(k+1,dim-1);
447 const unsigned numInnerEdgeDofs = k-1;
468 subEntity[lastIndex++] = 0;
469 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
470 subEntity[lastIndex++] = 6;
472 subEntity[lastIndex++] = 1;
475 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
476 subEntity[lastIndex++] = 4;
477 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
478 subEntity[lastIndex++] = 4;
479 subEntity[lastIndex++] = 5;
483 subEntity[lastIndex++] = 2;
484 for (
unsigned i = 0; i < k - 1; ++i)
485 subEntity[lastIndex++] = 7;
486 subEntity[lastIndex++] = 3;
488 assert(numFaceIndices==lastIndex);
492 for(
unsigned f = 0; f < numInnerEdgeDofs; ++f) {
495 subEntity[lastIndex++] = 0;
496 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
497 subEntity[lastIndex++] = 2;
498 subEntity[lastIndex++] = 1;
501 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
502 subEntity[lastIndex++] = 0;
503 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
504 subEntity[lastIndex++] = 0;
505 subEntity[lastIndex++] = 1;
509 subEntity[lastIndex++] = 2;
510 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
511 subEntity[lastIndex++] = 3;
512 subEntity[lastIndex++] = 3;
514 assert(lastIndex==(f+1+1)*numFaceIndices);
519 subEntity[lastIndex++] = 4;
520 for (
unsigned i = 0; i < k - 1; ++i)
521 subEntity[lastIndex++] = 10;
522 subEntity[lastIndex++] = 5;
525 for (
unsigned e = 0; e < k - 1; ++e) {
526 subEntity[lastIndex++] = 8;
527 for (
unsigned i = 0; i < k - 1; ++i)
528 subEntity[lastIndex++] = 5;
529 subEntity[lastIndex++] = 9;
533 subEntity[lastIndex++] = 6;
534 for (
unsigned i = 0; i < k - 1; ++i)
535 subEntity[lastIndex++] = 11;
536 subEntity[lastIndex++] = 7;
538 assert(numIndices==lastIndex);
543 LagrangeCubeLocalCoefficients ()
548 localKeys_[0] = LocalKey(0,0,0);
554 for (std::size_t i=0; i<size(); i++)
555 localKeys_[i] = LocalKey(i,dim,0);
562 std::vector<unsigned int> codim(size());
564 for (std::size_t i=0; i<codim.size(); i++)
570 std::array<unsigned int,dim> mIdx = multiindex(i);
571 for (
unsigned int j=0; j<dim; j++)
572 if (mIdx[j]==0 or mIdx[j]==k)
581 std::vector<unsigned int> index(size());
583 for (std::size_t i=0; i<size(); i++)
587 std::array<unsigned int,dim> mIdx = multiindex(i);
589 for (
int j=dim-1; j>=0; j--)
590 if (mIdx[j]>0 && mIdx[j]<k)
591 index[i] = (k-1)*index[i] + (mIdx[j]-1);
595 std::vector<unsigned int> subEntity(size());
612 for (
size_t i=0; i<size(); i++)
613 localKeys_[i] = LocalKey(subEntity[i], codim[i], index[i]);
617 static constexpr std::size_t size ()
619 return power(k+1,dim);
623 const LocalKey& localKey (std::size_t i)
const
625 return localKeys_[i];
629 std::vector<LocalKey> localKeys_;
636 template<
class LocalBasis>
637 class LagrangeCubeLocalInterpolation
648 template<
typename F,
typename C>
649 void interpolate (
const F& ff, std::vector<C>& out)
const
651 constexpr auto dim = LocalBasis::Traits::dimDomain;
652 constexpr auto k = LocalBasis::order();
653 using D =
typename LocalBasis::Traits::DomainFieldType;
655 typename LocalBasis::Traits::DomainType x;
656 auto&& f = Impl::makeFunctionWithCallOperator<typename LocalBasis::Traits::DomainType>(ff);
658 out.resize(LocalBasis::size());
671 for (
unsigned int i=0; i<LocalBasis::size(); i++)
674 for (
int j=0; j<dim; j++)
675 x[j] = (i & (1<<j)) ? 1.0 : 0.0;
683 for (
unsigned int i=0; i<LocalBasis::size(); i++)
686 std::array<unsigned int,dim> alpha(LocalBasis::multiindex(i));
689 for (
unsigned int j=0; j<dim; j++)
690 x[j] = (1.0*alpha[j])/k;
709 template<
class D,
class R,
int dim,
int k>
716 Impl::LagrangeCubeLocalCoefficients<dim,k>,
717 Impl::LagrangeCubeLocalInterpolation<Impl::LagrangeCubeLocalBasis<D,R,dim,k> > >;
737 return coefficients_;
744 return interpolation_;
748 static constexpr std::size_t
size ()
750 return power(k+1,dim);
761 Impl::LagrangeCubeLocalBasis<D,R,dim,k> basis_;
762 Impl::LagrangeCubeLocalCoefficients<dim,k> coefficients_;
763 Impl::LagrangeCubeLocalInterpolation<Impl::LagrangeCubeLocalBasis<D,R,dim,k> > interpolation_;
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:125
Lagrange finite element for cubes with arbitrary compile-time dimension and polynomial order.
Definition: lagrangecube.hh:711
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: lagrangecube.hh:728
LagrangeCubeLocalFiniteElement()
Default constructor.
Definition: lagrangecube.hh:724
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: lagrangecube.hh:742
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: lagrangecube.hh:755
static constexpr std::size_t size()
The number of shape functions.
Definition: lagrangecube.hh:748
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: lagrangecube.hh:735
Default exception for dummy implementations.
Definition: exceptions.hh:263
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.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:472
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:291
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
static const ReferenceElement & cube()
get hypercube reference elements
Definition: referenceelements.hh:210
D DomainType
domain type
Definition: localbasis.hh:42
traits helper struct
Definition: localfiniteelementtraits.hh:13
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24