5#ifndef DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGESIMPLEX_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGESIMPLEX_HH
19#include <dune/common/hybridutilities.hh>
22#include <dune/common/std/mdarray.hh>
23#include <dune/common/std/mdspan.hh>
25#include <dune/geometry/referenceelements.hh>
27#include <dune/localfunctions/common/localbasis.hh>
28#include <dune/localfunctions/common/localfiniteelementtraits.hh>
29#include <dune/localfunctions/common/localkey.hh>
31namespace Dune {
namespace Impl
34 template<
unsigned int dim,
int compileTimeOrder>
35 struct LagrangeSimplexTraits;
38 template<
unsigned int dim,
int compileTimeOrder>
39 requires (compileTimeOrder >= 0)
40 struct LagrangeSimplexTraits<dim,compileTimeOrder>
42 static constexpr bool is_static_order =
true;
44 constexpr LagrangeSimplexTraits(
int = compileTimeOrder) {}
47 static constexpr unsigned int order() {
return compileTimeOrder; }
53 static constexpr unsigned int size() {
return binomial(compileTimeOrder+dim,dim); }
58 template<
unsigned int dim>
59 struct LagrangeSimplexTraits<dim,-1>
64 static constexpr bool is_static_order =
false;
66 constexpr explicit LagrangeSimplexTraits(
int runTimeOrder)
67 : order_(runTimeOrder >= 0 ? (unsigned int)(runTimeOrder) : 0u)
75 constexpr unsigned int order()
const
81 constexpr std::size_t
size()
const
96 template<
class R,
unsigned int dim,
int k>
97 struct LagrangeSimplexLocalBasisBuffers
99 LagrangeSimplexLocalBasisBuffers(
int ) {}
102 static auto makeValueBuffer(
int )
104 using E = Std::extents<int,dim+1,k+1>;
105 using C = std::array<R,(dim+1)*(k+1)>;
106 return Std::mdarray<R,E,Std::layout_right,C>{};
110 static auto makeJacobianBuffer(
int )
112 using E = Std::extents<int,dim+1,2,k+1>;
113 using C = std::array<R,(dim+1)*2*(k+1)>;
114 return Std::mdarray<R,E,Std::layout_right,C>{};
118 template <
class T, T partialOrders>
119 static auto makePartialBuffer(std::integral_constant<T,partialOrders>,
int )
121 using E = Std::extents<int,dim+1,partialOrders+1,k+1>;
122 using C = std::array<R,(dim+1)*(partialOrders+1)*(k+1)>;
123 return Std::mdarray<R,E,Std::layout_right,C>{};
135 template<
class R,
unsigned int dim>
136 struct LagrangeSimplexLocalBasisBuffers<R,dim,-1>
138 LagrangeSimplexLocalBasisBuffers(
int order)
139 : buffer_((dim+1)*
std::
max<int>(2,order+1)*(order+1))
142 mutable std::vector<R> buffer_{};
145 auto makeValueBuffer(
int order)
const
147 using E = Std::extents<int,dim+1,std::dynamic_extent>;
148 return Std::mdspan<R,E>{buffer_.data(), E(order+1)};
152 auto makeJacobianBuffer(
int order)
const
154 using E = Std::extents<int,dim+1,2,std::dynamic_extent>;
155 return Std::mdspan<R,E>{buffer_.data(), E(order+1)};
159 auto makePartialBuffer(
int partialOrders,
int order)
const
161 using E = Std::extents<int,dim+1,std::dynamic_extent,std::dynamic_extent>;
162 return Std::mdspan<R,E>{buffer_.data(), E(partialOrders+1, order+1)};
177 template<
class D,
class R,
unsigned int dim,
int compileTimeOrder>
178 class LagrangeSimplexLocalBasis
179 :
private LagrangeSimplexTraits<dim,compileTimeOrder>
180 ,
private LagrangeSimplexLocalBasisBuffers<R,dim,compileTimeOrder>
182 using OrderTraits = LagrangeSimplexTraits<dim,compileTimeOrder>;
183 using Buffers = LagrangeSimplexLocalBasisBuffers<R,dim,compileTimeOrder>;
184 static constexpr bool is_static_order = OrderTraits::is_static_order;
188 constexpr LagrangeSimplexLocalBasis ()
requires (OrderTraits::is_static_order)
189 : LagrangeSimplexLocalBasis(OrderTraits())
193 explicit constexpr LagrangeSimplexLocalBasis(OrderTraits orderTraits)
194 : OrderTraits(orderTraits)
195 , Buffers(orderTraits.order())
198 using OrderTraits::order;
199 using OrderTraits::size;
206 template <
class I, std::size_t e0, std::size_t... ee,
class A>
207 static auto slice(Std::mdspan<R,Std::extents<I,e0,ee...>,Std::layout_right,A> L,
int i)
210 return Std::mdspan<R,Std::extents<I,ee...>,Std::layout_right>{
211 L.data_handle() + i * L.stride(0), L.extent(ii+1)...};
212 }, std::make_index_sequence<
sizeof...(ee)>{});
216 template <
class I, std::size_t e0, std::size_t... ee,
class C>
217 static auto slice(Std::mdarray<R,Std::extents<I,e0,ee...>,Std::layout_right,C>& L,
int i)
219 return slice(L.to_mdspan(), i);
229 constexpr auto barycentric(
const auto& x)
const
231 auto b = std::array<R,dim+1>{};
235 b[i] = order() * x[i];
245 static constexpr void evaluateLagrangePolynomials(
const R& t,
auto&& L,
int k)
249 L[i+1] = L[i] * (t - i) / (i+1);
254 static constexpr void evaluateLagrangePolynomialDerivative(
const R& t,
auto&& LL,
unsigned int maxDerivativeOrder,
int k)
256 auto L = slice(LL,0);
259 L[i+1] = L[i] * (t - i) / (i+1);
262 auto F = slice(LL,j);
263 auto DF = slice(LL,j+1);
266 DF[i+1] = (DF[i] * (t - i) + (j+1)*F[i]) / (i+1);
270 using BarycentricMultiIndex = std::array<int,dim+1>;
286 static constexpr R barycentricDerivative(
287 BarycentricMultiIndex beta,
290 const BarycentricMultiIndex& i,
291 const BarycentricMultiIndex& alpha = {})
304 auto leftDerivatives = alpha;
305 auto rightDerivatives = alpha;
306 leftDerivatives[j]++;
307 rightDerivatives.back()++;
309 return (barycentricDerivative(beta, L, k, i, leftDerivatives) -
310 barycentricDerivative(beta, L, k, i, rightDerivatives)) * k;
320 y *= L(j,alpha[j],i[j]);
326 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
330 std::vector<typename Traits::RangeType>& out)
const
332 const int k = order();
346 for (
size_t i=0; i<dim; i++)
355 auto z = barycentric(x);
357 auto L = Buffers::makeValueBuffer(k);
359 evaluateLagrangePolynomials(z[j], slice(L,j), k);
365 for (
auto i1 : std::array{k - i0})
366 out[n++] = L(0,i0) * L(1,i1);
374 for (
auto i2 : std::array{k - i1 - i0})
375 out[n++] = L(0,i0) * L(1,i1) * L(2,i2);
384 for (
auto i3 : std::array{k - i2 - i1 - i0})
385 out[n++] = L(0,i0) * L(1,i1) * L(2,i2) * L(3,i3);
389 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
398 std::vector<typename Traits::JacobianType>& out)
const
400 const int k = order();
406 std::fill(out[0][0].begin(), out[0][0].end(), 0);
413 std::fill(out[0][0].begin(), out[0][0].end(), -1);
415 for (
unsigned int i=0; i<dim; i++)
416 for (
unsigned int j=0; j<dim; j++)
417 out[i+1][0][j] = (i==j);
423 auto z = barycentric(x);
426 auto L = Buffers::makeJacobianBuffer(k);
428 evaluateLagrangePolynomialDerivative(z[j], slice(L,j), 1, k);
435 for (
auto i1 : std::array{k-i0})
437 out[n][0][0] = (L(0,1,i0) * L(1,0,i1) - L(0,0,i0) * L(1,1,i1))*k;
450 for (
auto i2 : std::array{k - i1 - i0})
452 out[n][0][0] = (L(0,1,i0) * L(1,0,i1) * L(2,0,i2) - L(0,0,i0) * L(1,0,i1) * L(2,1,i2))*k;
453 out[n][0][1] = (L(0,0,i0) * L(1,1,i1) * L(2,0,i2) - L(0,0,i0) * L(1,0,i1) * L(2,1,i2))*k;
469 for (
auto i3 : std::array{k - i2 - i1 - i0})
471 out[n][0][0] = (L(0,1,i0) * L(1,0,i1) * L(2,0,i2) * L(3,0,i3) - L(0,0,i0) * L(1,0,i1) * L(2,0,i2) * L(3,1,i3))*k;
472 out[n][0][1] = (L(0,0,i0) * L(1,1,i1) * L(2,0,i2) * L(3,0,i3) - L(0,0,i0) * L(1,0,i1) * L(2,0,i2) * L(3,1,i3))*k;
473 out[n][0][2] = (L(0,0,i0) * L(1,0,i1) * L(2,1,i2) * L(3,0,i3) - L(0,0,i0) * L(1,0,i1) * L(2,0,i2) * L(3,1,i3))*k;
483 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
492 void partial(
const std::array<unsigned int,dim>& partialOrders,
494 std::vector<typename Traits::RangeType>& out)
const
496 const int k = order();
497 int totalOrder =
std::accumulate(partialOrders.begin(), partialOrders.end(), 0u);
504 evaluateFunction(in, out);
511 for(
auto& out_i : out)
522 auto direction = std::find(partialOrders.begin(), partialOrders.end(), 1);
524 for (
unsigned int i=0; i<dim; i++)
525 out[i+1] = (i==(direction-partialOrders.begin()));
533 auto supportedOrders = [&]{
534 if constexpr(is_static_order)
542 auto z = barycentric(in);
545 auto L = Buffers::makePartialBuffer(staticTotalOrder,k);
547 evaluateLagrangePolynomialDerivative(z[j], slice(L,j), partialOrders[j], k);
548 evaluateLagrangePolynomialDerivative(z[dim], slice(L,dim), totalOrder, k);
550 auto barycentricOrder = BarycentricMultiIndex{};
552 barycentricOrder[j] = partialOrders[j];
553 barycentricOrder[dim] = 0;
555 if constexpr (dim==1)
559 for (
auto i1 : std::array{k - i0})
560 out[n++] = barycentricDerivative(barycentricOrder, L, k, BarycentricMultiIndex{i0, i1});
562 if constexpr (dim==2)
567 for (
auto i2 : std::array{k - i1 - i0})
568 out[n++] = barycentricDerivative(barycentricOrder, L, k, BarycentricMultiIndex{i0, i1, i2});
570 if constexpr (dim==3)
576 for (
auto i3 : std::array{k - i2 - i1 - i0})
577 out[n++] = barycentricDerivative(barycentricOrder, L, k, BarycentricMultiIndex{i0, i1, i2, i3});
588 template<
unsigned int dim,
int compileTimeOrder>
589 class LagrangeSimplexLocalCoefficients
590 :
private LagrangeSimplexTraits<dim,compileTimeOrder>
592 using OrderTraits = LagrangeSimplexTraits<dim,compileTimeOrder>;
596 LagrangeSimplexLocalCoefficients ()
requires (OrderTraits::is_static_order)
597 : LagrangeSimplexLocalCoefficients(OrderTraits())
601 LagrangeSimplexLocalCoefficients (OrderTraits orderTraits)
602 : OrderTraits(orderTraits)
603 , localKeys_(OrderTraits::
size())
605 const int k = OrderTraits::order();
608 localKeys_[0] = LocalKey(0,0,0);
614 for (std::size_t i=0; i<OrderTraits::size(); i++)
615 localKeys_[i] = LocalKey(i,dim,0);
622 localKeys_[0] = LocalKey(0,1,0);
623 for (
int i=1; i<k; i++)
624 localKeys_[i] = LocalKey(0,0,i-1);
625 localKeys_.back() = LocalKey(1,1,0);
633 for (
int j=0; j<=k; j++)
634 for (
int i=0; i<=k-j; i++)
638 localKeys_[n++] = LocalKey(0,2,0);
643 localKeys_[n++] = LocalKey(1,2,0);
648 localKeys_[n++] = LocalKey(2,2,0);
653 localKeys_[n++] = LocalKey(0,1,i-1);
658 localKeys_[n++] = LocalKey(1,1,j-1);
663 localKeys_[n++] = LocalKey(2,1,j-1);
666 localKeys_[n++] = LocalKey(0,0,c++);
673 std::array<unsigned int, dim+1> vertexMap;
674 for (
unsigned int i=0; i<=dim; i++)
676 generateLocalKeys(vertexMap);
679 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for k<=1 or dim<=3!");
688 constexpr LagrangeSimplexLocalCoefficients (
const std::array<unsigned int, dim+1> vertexMap)
requires (OrderTraits::is_static_order)
689 : localKeys_(OrderTraits::size())
691 if (dim!=2 && dim!=3)
692 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
694 generateLocalKeys(vertexMap);
698 template<std::input_iterator VertexMap>
699 constexpr LagrangeSimplexLocalCoefficients(
const VertexMap &vertexmap)
requires (OrderTraits::is_static_order)
700 : localKeys_(OrderTraits::size())
702 if (dim!=2 && dim!=3)
703 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
705 std::array<unsigned int, dim+1> vertexmap_array;
706 std::copy(vertexmap, vertexmap + dim + 1, vertexmap_array.begin());
707 generateLocalKeys(vertexmap_array);
711 constexpr const LocalKey& localKey (std::size_t i)
const
713 return localKeys_[i];
716 using OrderTraits::order;
717 using OrderTraits::size;
720 std::vector<LocalKey> localKeys_;
722 constexpr void generateLocalKeys(
const std::array<unsigned int, dim+1> vertexMap)
724 const int k = OrderTraits::order();
727 localKeys_[0] = LocalKey(0,0,0);
736 for (
int j=0; j<=k; j++)
737 for (
int i=0; i<=k-j; i++)
741 localKeys_[n++] = LocalKey(0,2,0);
746 localKeys_[n++] = LocalKey(1,2,0);
751 localKeys_[n++] = LocalKey(2,2,0);
756 localKeys_[n++] = LocalKey(0,1,i-1);
761 localKeys_[n++] = LocalKey(1,1,j-1);
766 localKeys_[n++] = LocalKey(2,1,j-1);
769 localKeys_[n++] = LocalKey(0,0,c++);
774 flip[0] = vertexMap[0] > vertexMap[1];
775 flip[1] = vertexMap[0] > vertexMap[2];
776 flip[2] = vertexMap[1] > vertexMap[2];
777 for (std::size_t i=0; i<OrderTraits::size(); i++)
778 if (localKeys_[i].codim()==1 && flip[localKeys_[i].subEntity()])
779 localKeys_[i].index(k-2-localKeys_[i].index());
785 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for dim==3!");
787 unsigned int subindex[16];
788 unsigned int codim_count[4] = {0};
789 for (
unsigned int m = 1; m < 16; ++m)
791 unsigned int codim = !(m&1) + !(m&2) + !(m&4) + !(m&8);
792 subindex[m] = codim_count[codim]++;
795 int a1 = (3*k + 12)*k + 11;
797 unsigned int dof_count[16] = {0};
799 for (i[3] = 0; i[3] <= k; ++i[3])
800 for (i[2] = 0; i[2] <= k - i[3]; ++i[2])
801 for (i[1] = 0; i[1] <= k - i[2] - i[3]; ++i[1])
803 i[0] = k - i[1] - i[2] - i[3];
805 unsigned int entity = 0;
806 unsigned int codim = 0;
807 for (
unsigned int m = 0; m < 4; ++m)
809 j[m] = i[vertexMap[m]];
810 entity += !!j[m] << m;
813 int local_index = j[3]*(a1 + (a2 + j[3])*j[3])/6
814 + j[2]*(2*(k - j[3]) + 3 - j[2])/2 + j[1];
815 localKeys_[local_index] = LocalKey(subindex[entity], codim, dof_count[entity]++);
827 template<
class D,
class R,
int dim,
int compileTimeOrder>
828 class LagrangeSimplexLocalInterpolation
829 :
private LagrangeSimplexTraits<dim, compileTimeOrder>
831 using OrderTraits = LagrangeSimplexTraits<dim, compileTimeOrder>;
832 using Traits =
typename LagrangeSimplexLocalBasis<D, R, dim, compileTimeOrder>::Traits;
834 using OrderTraits::size;
837 constexpr LagrangeSimplexLocalInterpolation ()
requires (OrderTraits::is_static_order)
842 explicit constexpr LagrangeSimplexLocalInterpolation (OrderTraits orderTraits)
843 : OrderTraits(orderTraits)
853 template<
typename F,
typename C>
854 constexpr void interpolate (
const F& f, std::vector<C>& out)
const
856 const int k = OrderTraits::order();
857 out.resize(OrderTraits::size());
867 typename Traits::DomainType x;
873 std::fill(x.begin(), x.end(), 0);
877 for (
unsigned int i=0; i<dim; i++)
879 for (
unsigned int j=0; j<dim; j++)
889 for (
int i=0; i<k+1; i++)
900 for (
int j=0; j<=k; j++)
901 for (
int i=0; i<=k-j; i++)
903 x = { ((D)i)/k, ((D)j)/k };
911 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalInterpolation only implemented for dim<=3!");
913 const int kdiv = (k==0 ? 1 : k);
916 for (
int i2 = 0; i2 <= k; i2++)
917 for (
int i1 = 0; i1 <= k-i2; i1++)
918 for (
int i0 = 0; i0 <= k-i1-i2; i0++)
920 x[0] = ((D)i0)/((D)kdiv);
921 x[1] = ((D)i1)/((D)kdiv);
922 x[2] = ((D)i2)/((D)kdiv);
987 template<
class D,
class R,
int dim,
int compileTimeOrder = -1>
989 :
private Impl::LagrangeSimplexTraits<dim, compileTimeOrder>
991 using OrderTraits = Impl::LagrangeSimplexTraits<dim, compileTimeOrder>;
997 Impl::LagrangeSimplexLocalBasis<D,R,dim,compileTimeOrder>,
998 Impl::LagrangeSimplexLocalCoefficients<dim,compileTimeOrder>,
999 Impl::LagrangeSimplexLocalInterpolation<D,R,dim,compileTimeOrder> >;
1005 , coefficients_(*this)
1006 , interpolation_(*this)
1008 static_assert(OrderTraits::is_static_order,
"Default constructor only allowed for compile-time order >= 0");
1009 static_assert(compileTimeOrder >= 0,
"Default constructor only allowed for compile-time order >= 0");
1014 : OrderTraits(runTimeOrder)
1016 , coefficients_(*this)
1017 , interpolation_(*this)
1019 if (runTimeOrder < 0)
1021 if (compileTimeOrder >= 0 && compileTimeOrder != runTimeOrder)
1027 template<
typename VertexMap>
1031 , coefficients_(vertexmap)
1032 , interpolation_(*
this)
1046 return coefficients_;
1053 return interpolation_;
1057 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 simplices with arbitrary compile-time dimension and static or dynamic pol...
Definition: lagrangesimplex.hh:990
constexpr const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition: lagrangesimplex.hh:1044
constexpr LagrangeSimplexLocalFiniteElement(const VertexMap &vertexmap)
Constructs a finite element given a vertex reordering.
Definition: lagrangesimplex.hh:1028
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition: lagrangesimplex.hh:1061
constexpr const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition: lagrangesimplex.hh:1037
constexpr const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition: lagrangesimplex.hh:1051
constexpr LagrangeSimplexLocalFiniteElement(int runTimeOrder)
Constructor for run-time order.
Definition: lagrangesimplex.hh:1013
constexpr LagrangeSimplexLocalFiniteElement()
Constructor for compile-time order.
Definition: lagrangesimplex.hh:1002
A few common exception classes.
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.
decltype(auto) constexpr unpackIntegerSequence(F &&f, std::integer_sequence< I, i... > sequence)
Unpack an std::integer_sequence<I,i...> to std::integral_constant<I,i>...
Definition: indices.hh:124
std::integral_constant< std::size_t, i > index_constant
An index constant with value i.
Definition: indices.hh:29
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:453
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr decltype(auto) switchCases(const Cases &cases, const Value &value, Branches &&branches, ElseBranch &&elseBranch)
Switch statement.
Definition: hybridutilities.hh:661
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:489
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:284
static constexpr IntegralRange< std::decay_t< T > > range(T &&from, U &&to) noexcept
free standing function for setting up a range based for loop over an integer range for (auto i: range...
Definition: rangeutilities.hh:288
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
static constexpr T binomial(const T &n, const T &k) noexcept
calculate the binomial coefficient n over k as a constexpr
Definition: math.hh:113
Utilities for reduction like operations on ranges.
static const ReferenceElement & simplex()
get simplex reference elements
Definition: referenceelements.hh:162
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