5#ifndef DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGECUBE_HH
16#include <dune/geometry/referenceelements.hh>
18#include <dune/localfunctions/common/localbasis.hh>
19#include <dune/localfunctions/common/localfiniteelementtraits.hh>
20#include <dune/localfunctions/common/localkey.hh>
22namespace Dune {
namespace Impl
25 template<
unsigned int dim,
int compileTimeOrder>
26 struct LagrangeCubeOrderTraits;
29 template<
unsigned int dim,
int compileTimeOrder>
30 requires (compileTimeOrder >= 0)
31 struct LagrangeCubeOrderTraits<dim,compileTimeOrder>
33 static constexpr bool is_static_order =
true;
35 constexpr LagrangeCubeOrderTraits(
int = compileTimeOrder) {}
40 static constexpr std::size_t
size ()
42 return power(compileTimeOrder+1, dim);
48 static constexpr unsigned int order ()
50 return compileTimeOrder;
54 static constexpr std::array<unsigned int,dim> multiindex (
unsigned int i)
56 std::array<unsigned int,dim> alpha;
57 for (
unsigned int j=0; j<dim; j++)
59 alpha[j] = i % (order()+1);
67 template<
unsigned int dim>
68 struct LagrangeCubeOrderTraits<dim,-1>
73 static constexpr bool is_static_order =
false;
78 constexpr LagrangeCubeOrderTraits (
int runTimeOrder)
79 : order_(runTimeOrder >= 0 ? (unsigned int)(runTimeOrder) : 0u)
80 , size_(
power(order_+1, dim))
89 constexpr std::size_t
size ()
const
97 constexpr unsigned int order ()
const
103 constexpr std::array<unsigned int,dim> multiindex (
unsigned int i)
const
105 std::array<unsigned int,dim> alpha;
106 for (
unsigned int j=0; j<dim; j++)
108 alpha[j] = i % (order()+1);
125 template<
class D,
class R,
unsigned int dim,
int compileTimeOrder>
126 class LagrangeCubeLocalBasis
127 :
private LagrangeCubeOrderTraits<dim,compileTimeOrder>
129 using OrderTraits = LagrangeCubeOrderTraits<dim,compileTimeOrder>;
132 constexpr R p (
unsigned int i, D x)
const
135 for (
unsigned int j=0; j<=order(); j++)
136 if (j!=i) result *= (order()*x-j)/((
int)i-(int)j);
141 constexpr R dp (
unsigned int i, D x)
const
145 for (
unsigned int j=0; j<=order(); j++)
149 R prod( (order()*1.0)/((
int)i-(
int)j) );
150 for (
unsigned int l=0; l<=order(); l++)
152 prod *= (order()*x-l)/((
int)i-(int)l);
161 constexpr R ddp(
unsigned int j, D x)
const
165 for (
unsigned int i=0; i<=order(); i++)
172 for (
unsigned int m=0; m<=order(); m++)
177 R prod( (order()*1.0)/((
int)j-(
int)m) );
178 for (
unsigned int l=0; l<=order(); l++)
179 if (l!=i && l!=j && l!=m)
180 prod *= (order()*x-l)/((
int)j-(int)l);
184 result += sum * ( (order()*1.0)/((
int)j-(int)i) );
190 using OrderTraits::multiindex;
193 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
195 constexpr LagrangeCubeLocalBasis ()
requires (OrderTraits::is_static_order)
201 explicit constexpr LagrangeCubeLocalBasis (OrderTraits orderTraits)
202 : OrderTraits(orderTraits)
205 using OrderTraits::size;
206 using OrderTraits::order;
210 std::vector<typename Traits::RangeType>& out)
const
223 for (
size_t i=0; i<
size(); i++)
227 for (
unsigned int j=0; j<dim; j++)
229 out[i] *= (i & (1<<j)) ? x[j] : 1-x[j];
235 for (
size_t i=0; i<
size(); i++)
238 std::array<unsigned int,dim> alpha(multiindex(i));
244 for (
unsigned int j=0; j<dim; j++)
245 out[i] *= p(alpha[j],x[j]);
255 std::vector<typename Traits::JacobianType>& out)
const
262 std::fill(out[0][0].begin(), out[0][0].end(), 0);
270 for (
size_t i=0; i<
size(); i++)
273 for (
unsigned int j=0; j<dim; j++)
277 out[i][0][j] = (i & (1<<j)) ? 1 : -1;
279 for (
unsigned int l=0; l<dim; l++)
283 out[i][0][j] *= (i & (1<<l)) ? x[l] : 1-x[l];
293 for (
size_t i=0; i<
size(); i++)
296 std::array<unsigned int,dim> alpha(multiindex(i));
299 for (
unsigned int j=0; j<dim; j++)
303 out[i][0][j] = dp(alpha[j],x[j]);
306 for (
unsigned int l=0; l<dim; l++)
308 out[i][0][j] *= p(alpha[l],x[l]);
319 constexpr void partial (
const std::array<unsigned int,dim>& partialOrders,
321 std::vector<typename Traits::RangeType>& out)
const
323 auto totalOrder =
std::accumulate(partialOrders.begin(), partialOrders.end(), 0);
329 out[0] = (totalOrder==0);
337 evaluateFunction(in, out);
339 else if (totalOrder == 1)
343 auto direction = std::distance(partialOrders.begin(),
344 std::find(partialOrders.begin(), partialOrders.end(), 1));
345 if (direction >= dim)
346 DUNE_THROW(RangeError,
"Direction of partial derivative not found!");
349 for (std::size_t i = 0; i <
size(); ++i)
353 out[i] = (i & (1<<direction)) ? 1 : -1;
355 for (
unsigned int j = 0; j < dim; ++j)
359 out[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
363 else if (totalOrder == 2)
366 for (
size_t i=0; i<
size(); i++)
369 std::array<unsigned int,dim> alpha(multiindex(i));
375 for (std::size_t l=0; l<dim; l++)
377 switch (partialOrders[l])
380 out[i][0] *= p(alpha[l],in[l]);
384 out[i][0] *= dp(alpha[l],in[l]);
390 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
396 DUNE_THROW(NotImplemented,
"Partial derivative of order " << totalOrder <<
" is not implemented!");
404 for (
size_t i=0; i<
size(); i++)
407 std::array<unsigned int,dim> alpha(multiindex(i));
413 for (std::size_t l=0; l<dim; l++)
415 switch (partialOrders[l])
418 out[i][0] *= p(alpha[l],in[l]);
421 out[i][0] *= dp(alpha[l],in[l]);
424 out[i][0] *= ddp(alpha[l],in[l]);
427 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
439 template<
unsigned int dim,
int compileTimeOrder>
440 class LagrangeCubeLocalCoefficients
441 :
private LagrangeCubeOrderTraits<dim,compileTimeOrder>
443 using OrderTraits = LagrangeCubeOrderTraits<dim,compileTimeOrder>;
445 using OrderTraits::multiindex;
448 void setup1d (std::vector<unsigned int>& subEntity)
450 const unsigned int k = order();
453 unsigned lastIndex=0;
460 subEntity[lastIndex++] = 0;
461 for (
unsigned i = 0; i < k - 1; ++i)
462 subEntity[lastIndex++] = 0;
464 subEntity[lastIndex++] = 1;
466 assert(
power(k+1,dim)==lastIndex);
469 void setup2d (std::vector<unsigned int>& subEntity)
471 const unsigned int k = order();
474 unsigned lastIndex=0;
488 subEntity[lastIndex++] = 0;
489 for (
unsigned i = 0; i < k - 1; ++i)
490 subEntity[lastIndex++] = 2;
492 subEntity[lastIndex++] = 1;
495 for (
unsigned e = 0; e < k - 1; ++e) {
496 subEntity[lastIndex++] = 0;
497 for (
unsigned i = 0; i < k - 1; ++i)
498 subEntity[lastIndex++] = 0;
499 subEntity[lastIndex++] = 1;
503 subEntity[lastIndex++] = 2;
504 for (
unsigned i = 0; i < k - 1; ++i)
505 subEntity[lastIndex++] = 3;
507 subEntity[lastIndex++] = 3;
509 assert(
power(k+1,dim)==lastIndex);
512 void setup3d (std::vector<unsigned int>& subEntity)
514 const unsigned int k = order();
517 unsigned lastIndex=0;
519 const unsigned numIndices =
power(k+1,dim);
520 const unsigned numFaceIndices =
power(k+1,dim-1);
522 const unsigned numInnerEdgeDofs = k-1;
543 subEntity[lastIndex++] = 0;
544 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
545 subEntity[lastIndex++] = 6;
547 subEntity[lastIndex++] = 1;
550 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
551 subEntity[lastIndex++] = 4;
552 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
553 subEntity[lastIndex++] = 4;
554 subEntity[lastIndex++] = 5;
558 subEntity[lastIndex++] = 2;
559 for (
unsigned i = 0; i < k - 1; ++i)
560 subEntity[lastIndex++] = 7;
561 subEntity[lastIndex++] = 3;
563 assert(numFaceIndices==lastIndex);
567 for(
unsigned f = 0; f < numInnerEdgeDofs; ++f) {
570 subEntity[lastIndex++] = 0;
571 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
572 subEntity[lastIndex++] = 2;
573 subEntity[lastIndex++] = 1;
576 for (
unsigned e = 0; e < numInnerEdgeDofs; ++e) {
577 subEntity[lastIndex++] = 0;
578 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
579 subEntity[lastIndex++] = 0;
580 subEntity[lastIndex++] = 1;
584 subEntity[lastIndex++] = 2;
585 for (
unsigned i = 0; i < numInnerEdgeDofs; ++i)
586 subEntity[lastIndex++] = 3;
587 subEntity[lastIndex++] = 3;
589 assert(lastIndex==(f+1+1)*numFaceIndices);
594 subEntity[lastIndex++] = 4;
595 for (
unsigned i = 0; i < k - 1; ++i)
596 subEntity[lastIndex++] = 10;
597 subEntity[lastIndex++] = 5;
600 for (
unsigned e = 0; e < k - 1; ++e) {
601 subEntity[lastIndex++] = 8;
602 for (
unsigned i = 0; i < k - 1; ++i)
603 subEntity[lastIndex++] = 5;
604 subEntity[lastIndex++] = 9;
608 subEntity[lastIndex++] = 6;
609 for (
unsigned i = 0; i < k - 1; ++i)
610 subEntity[lastIndex++] = 11;
611 subEntity[lastIndex++] = 7;
613 assert(numIndices==lastIndex);
620 const unsigned int k = order();
621 localKeys_.resize(
size());
625 localKeys_[0] = LocalKey(0,0,0);
631 for (std::size_t i=0; i<
size(); i++)
632 localKeys_[i] = LocalKey(i,dim,0);
639 std::vector<unsigned int> codim(
size());
641 for (std::size_t i=0; i<codim.size(); i++)
647 std::array<unsigned int,dim> mIdx = multiindex(i);
648 for (
unsigned int j=0; j<dim; j++)
649 if (mIdx[j]==0 or mIdx[j]==k)
658 std::vector<unsigned int> index(
size());
660 for (std::size_t i=0; i<
size(); i++)
664 std::array<unsigned int,dim> mIdx = multiindex(i);
666 for (
int j=dim-1; j>=0; j--)
667 if (mIdx[j]>0 && mIdx[j]<k)
668 index[i] = (k-1)*index[i] + (mIdx[j]-1);
672 std::vector<unsigned int> subEntity(
size());
689 for (
size_t i=0; i<
size(); i++)
690 localKeys_[i] = LocalKey(subEntity[i], codim[i], index[i]);
695 explicit LagrangeCubeLocalCoefficients (OrderTraits orderTraits)
696 : OrderTraits(orderTraits)
701 using OrderTraits::size;
702 using OrderTraits::order;
705 constexpr const LocalKey& localKey (std::size_t i)
const
707 return localKeys_[i];
711 std::vector<LocalKey> localKeys_;
721 template<
class D,
class R,
unsigned int dim,
int compileTimeOrder>
722 class LagrangeCubeLocalInterpolation
723 :
private LagrangeCubeOrderTraits<dim, compileTimeOrder>
725 using OrderTraits = LagrangeCubeOrderTraits<dim, compileTimeOrder>;
726 using Traits =
typename LagrangeCubeLocalBasis<D,R,dim,compileTimeOrder>::Traits;
728 using OrderTraits::multiindex;
729 using OrderTraits::size;
730 using OrderTraits::order;
734 explicit constexpr LagrangeCubeLocalInterpolation (OrderTraits orderTraits)
735 : OrderTraits(orderTraits)
745 template<
typename F,
typename C>
746 constexpr void interpolate (
const F& f, std::vector<C>& out)
const
748 const unsigned int k = order();
759 typename Traits::DomainType x;
764 for (
unsigned int i=0; i<
size(); i++)
767 for (
unsigned int j=0; j<dim; j++)
768 x[j] = (i & (1<<j)) ? 1.0 : 0.0;
776 for (
unsigned int i=0; i<
size(); i++)
779 std::array<unsigned int,dim> alpha(multiindex(i));
782 for (
unsigned int j=0; j<dim; j++)
783 x[j] = (1.0*alpha[j])/k;
802 template<
class D,
class R,
int dim,
int compileTimeOrder = -1>
804 :
private Impl::LagrangeCubeOrderTraits<dim, compileTimeOrder>
806 using OrderTraits = Impl::LagrangeCubeOrderTraits<dim, compileTimeOrder>;
812 Impl::LagrangeCubeLocalBasis<D,R,dim,compileTimeOrder>,
813 Impl::LagrangeCubeLocalCoefficients<dim,compileTimeOrder>,
814 Impl::LagrangeCubeLocalInterpolation<D,R,dim,compileTimeOrder> >;
820 , coefficients_(*this)
821 , interpolation_(*this)
823 static_assert(OrderTraits::is_static_order,
"Default constructor only allowed for compile-time order >= 0");
824 static_assert(compileTimeOrder >= 0,
"Default constructor only allowed for compile-time order >= 0");
829 : OrderTraits(runTimeOrder)
831 , coefficients_(*this)
832 , interpolation_(*this)
834 if (runTimeOrder < 0)
836 if (compileTimeOrder >= 0 && compileTimeOrder != runTimeOrder)
851 return coefficients_;
858 return interpolation_;
862 using OrderTraits::size;
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:375
Lagrange finite element for cubes with arbitrary compile-time dimension and polynomial order.
Definition: lagrangecube.hh:805
constexpr const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: lagrangecube.hh:856
constexpr LagrangeCubeLocalFiniteElement(int runTimeOrder)
Constructor for dynamic order.
Definition: lagrangecube.hh:828
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: lagrangecube.hh:866
constexpr const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: lagrangecube.hh:842
constexpr const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: lagrangecube.hh:849
constexpr LagrangeCubeLocalFiniteElement()
Constructor for static order.
Definition: lagrangecube.hh:817
Default exception for dummy implementations.
Definition: exceptions.hh:357
Traits for type conversions and type information.
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,...)
Definition: exceptions.hh:314
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:462
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:284
Some useful basic math stuff.
Dune namespace
Definition: alignedallocator.hh:13
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
static const ReferenceElement & cube()
get hypercube reference elements
Definition: referenceelements.hh:168
D DomainType
domain type
Definition: localbasis.hh:43
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