5#ifndef DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGESIMPLEX_HH
6#define DUNE_LOCALFUNCTIONS_LAGRANGE_LAGRANGESIMPLEX_HH
12#include <dune/common/exceptions.hh>
13#include <dune/common/fmatrix.hh>
14#include <dune/common/fvector.hh>
15#include <dune/common/hybridutilities.hh>
16#include <dune/common/math.hh>
17#include <dune/common/rangeutilities.hh>
19#include <dune/geometry/referenceelements.hh>
25namespace Dune {
namespace Impl
37 template<
class D,
class R,
unsigned int dim,
unsigned int k>
38 class LagrangeSimplexLocalBasis
48 constexpr auto barycentric(
const auto& x)
const
50 auto b = std::array<R,dim+1>{};
52 for(
auto i : Dune::range(dim))
64 static constexpr void evaluateLagrangePolynomials(
const R& t,
auto& L)
67 for (
auto i : Dune::range(k))
68 L[i+1] = L[i] * (t - i) / (i+1);
73 static constexpr void evaluateLagrangePolynomialDerivative(
const R& t,
auto& LL,
unsigned int maxDerivativeOrder)
77 for (
auto i : Dune::range(k))
78 L[i+1] = L[i] * (t - i) / (i+1);
79 for(
auto j : Dune::range(maxDerivativeOrder))
84 for (
auto i : Dune::range(k))
85 DF[i+1] = (DF[i] * (t - i) + (j+1)*F[i]) / (i+1);
89 using BarycentricMultiIndex = std::array<unsigned int,dim+1>;
105 static constexpr R barycentricDerivative(
106 BarycentricMultiIndex beta,
108 const BarycentricMultiIndex& i,
109 const BarycentricMultiIndex& alpha = {})
118 for(
auto j : Dune::range(dim))
122 auto leftDerivatives = alpha;
123 auto rightDerivatives = alpha;
124 leftDerivatives[j]++;
125 rightDerivatives.back()++;
127 return (barycentricDerivative(beta, L, i, leftDerivatives) - barycentricDerivative(beta, L, i, rightDerivatives))*k;
136 for(
auto j : Dune::range(dim+1))
137 y *= L[j][alpha[j]][i[j]];
143 using Traits = LocalBasisTraits<D,dim,FieldVector<D,dim>,R,1,FieldVector<R,1>,FieldMatrix<R,1,dim> >;
149 static constexpr unsigned int size ()
151 return binomial(k+dim,dim);
156 std::vector<typename Traits::RangeType>& out)
const
171 for (
size_t i=0; i<dim; i++)
180 auto z = barycentric(x);
182 auto L = std::array<std::array<R,k+1>, dim+1>();
183 for (
auto j : Dune::range(dim+1))
184 evaluateLagrangePolynomials(z[j], L[j]);
189 for (
auto i0 : Dune::range(k + 1))
190 for (
auto i1 : std::array{k - i0})
191 out[n++] = L[0][i0] * L[1][i1];
197 for (
auto i1 : Dune::range(k + 1))
198 for (
auto i0 : Dune::range(k - i1 + 1))
199 for (
auto i2 : std::array{k - i1 - i0})
200 out[n++] = L[0][i0] * L[1][i1] * L[2][i2];
206 for (
auto i2 : Dune::range(k + 1))
207 for (
auto i1 : Dune::range(k - i2 + 1))
208 for (
auto i0 : Dune::range(k - i2 - i1 + 1))
209 for (
auto i3 : std::array{k - i2 - i1 - i0})
210 out[n++] = L[0][i0] * L[1][i1] * L[2][i2] * L[3][i3];
214 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
223 std::vector<typename Traits::JacobianType>& out)
const
230 std::fill(out[0][0].begin(), out[0][0].end(), 0);
237 std::fill(out[0][0].begin(), out[0][0].end(), -1);
239 for (
unsigned int i=0; i<dim; i++)
240 for (
unsigned int j=0; j<dim; j++)
241 out[i+1][0][j] = (i==j);
247 auto z = barycentric(x);
250 auto L = std::array<std::array<std::array<R,k+1>, 2>, dim+1>();
251 for (
auto j : Dune::range(dim+1))
252 evaluateLagrangePolynomialDerivative(z[j], L[j], 1);
257 for (
auto i0 : Dune::range(k + 1))
259 for (
auto i1 : std::array{k-i0})
261 out[n][0][0] = (L[0][1][i0] * L[1][0][i1] - L[0][0][i0] * L[1][1][i1])*k;
270 for (
auto i1 : Dune::range(k + 1))
272 for (
auto i0 : Dune::range(k - i1 + 1))
274 for (
auto i2 : std::array{k - i1 - i0})
276 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;
277 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;
287 for (
auto i2 : Dune::range(k + 1))
289 for (
auto i1 : Dune::range(k - i2 + 1))
291 for (
auto i0 : Dune::range(k - i2 - i1 + 1))
293 for (
auto i3 : std::array{k - i2 - i1 - i0})
295 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;
296 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;
297 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;
307 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalBasis for k>=2 only implemented for dim<=3");
316 void partial(
const std::array<unsigned int,dim>& order,
318 std::vector<typename Traits::RangeType>& out)
const
320 auto totalOrder = std::accumulate(order.begin(), order.end(), 0u);
327 evaluateFunction(in, out);
334 for(
auto& out_i : out)
345 auto direction = std::find(order.begin(), order.end(), 1);
347 for (
unsigned int i=0; i<dim; i++)
348 out[i+1] = (i==(direction-order.begin()));
356 auto supportedStaticOrders = Dune::range(Dune::index_constant<1>{}, Dune::index_constant<k+1>{});
357 return Dune::Hybrid::switchCases(supportedStaticOrders, totalOrder, [&](
auto staticTotalOrder) {
360 auto z = barycentric(in);
363 auto L = std::array<std::array<std::array<R, k+1>, staticTotalOrder+1>, dim+1>();
364 for (
auto j : Dune::range(dim))
365 evaluateLagrangePolynomialDerivative(z[j], L[j], order[j]);
366 evaluateLagrangePolynomialDerivative(z[dim], L[dim], totalOrder);
368 auto barycentricOrder = BarycentricMultiIndex{};
369 for (
auto j : Dune::range(dim))
370 barycentricOrder[j] = order[j];
371 barycentricOrder[dim] = 0;
373 if constexpr (dim==1)
376 for (
auto i0 : Dune::range(k + 1))
377 for (
auto i1 : std::array{k - i0})
378 out[n++] = barycentricDerivative(barycentricOrder, L, BarycentricMultiIndex{i0, i1});
380 if constexpr (dim==2)
383 for (
auto i1 : Dune::range(k + 1))
384 for (
auto i0 : Dune::range(k - i1 + 1))
385 for (
auto i2 : std::array{k - i1 - i0})
386 out[n++] = barycentricDerivative(barycentricOrder, L, BarycentricMultiIndex{i0, i1, i2});
388 if constexpr (dim==3)
391 for (
auto i2 : Dune::range(k + 1))
392 for (
auto i1 : Dune::range(k - i2 + 1))
393 for (
auto i0 : Dune::range(k - i2 - i1 + 1))
394 for (
auto i3 : std::array{k - i2 - i1 - i0})
395 out[n++] = barycentricDerivative(barycentricOrder, L, BarycentricMultiIndex{i0, i1, i2, i3});
401 static constexpr unsigned int order ()
412 template<
unsigned int dim,
unsigned int k>
413 class LagrangeSimplexLocalCoefficients
417 LagrangeSimplexLocalCoefficients ()
422 localKeys_[0] = LocalKey(0,0,0);
428 for (std::size_t i=0; i<size(); i++)
429 localKeys_[i] = LocalKey(i,dim,0);
436 localKeys_[0] = LocalKey(0,1,0);
437 for (
unsigned int i=1; i<k; i++)
438 localKeys_[i] = LocalKey(0,0,i-1);
439 localKeys_.back() = LocalKey(1,1,0);
447 for (
unsigned int j=0; j<=k; j++)
448 for (
unsigned int i=0; i<=k-j; i++)
452 localKeys_[n++] = LocalKey(0,2,0);
457 localKeys_[n++] = LocalKey(1,2,0);
462 localKeys_[n++] = LocalKey(2,2,0);
467 localKeys_[n++] = LocalKey(0,1,i-1);
472 localKeys_[n++] = LocalKey(1,1,j-1);
477 localKeys_[n++] = LocalKey(2,1,j-1);
480 localKeys_[n++] = LocalKey(0,0,c++);
487 std::array<unsigned int, dim+1> vertexMap;
488 for (
unsigned int i=0; i<=dim; i++)
490 generateLocalKeys(vertexMap);
493 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for k<=1 or dim<=3!");
502 LagrangeSimplexLocalCoefficients (
const std::array<unsigned int, dim+1> vertexMap)
505 if (dim!=2 && dim!=3)
506 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
508 generateLocalKeys(vertexMap);
512 template<
class VertexMap>
513 LagrangeSimplexLocalCoefficients(
const VertexMap &vertexmap)
516 if (dim!=2 && dim!=3)
517 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for dim==2 and dim==3!");
519 std::array<unsigned int, dim+1> vertexmap_array;
520 std::copy(vertexmap, vertexmap + dim + 1, vertexmap_array.begin());
521 generateLocalKeys(vertexmap_array);
525 static constexpr std::size_t size ()
527 return binomial(k+dim,dim);
531 const LocalKey& localKey (std::size_t i)
const
533 return localKeys_[i];
537 std::vector<LocalKey> localKeys_;
539 void generateLocalKeys(
const std::array<unsigned int, dim+1> vertexMap)
543 localKeys_[0] = LocalKey(0,0,0);
552 for (
unsigned int j=0; j<=k; j++)
553 for (
unsigned int i=0; i<=k-j; i++)
557 localKeys_[n++] = LocalKey(0,2,0);
562 localKeys_[n++] = LocalKey(1,2,0);
567 localKeys_[n++] = LocalKey(2,2,0);
572 localKeys_[n++] = LocalKey(0,1,i-1);
577 localKeys_[n++] = LocalKey(1,1,j-1);
582 localKeys_[n++] = LocalKey(2,1,j-1);
585 localKeys_[n++] = LocalKey(0,0,c++);
590 flip[0] = vertexMap[0] > vertexMap[1];
591 flip[1] = vertexMap[0] > vertexMap[2];
592 flip[2] = vertexMap[1] > vertexMap[2];
593 for (std::size_t i=0; i<size(); i++)
594 if (localKeys_[i].codim()==1 && flip[localKeys_[i].subEntity()])
595 localKeys_[i].index(k-2-localKeys_[i].index());
601 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalCoefficients only implemented for dim==3!");
603 unsigned int subindex[16];
604 unsigned int codim_count[4] = {0};
605 for (
unsigned int m = 1; m < 16; ++m)
607 unsigned int codim = !(m&1) + !(m&2) + !(m&4) + !(m&8);
608 subindex[m] = codim_count[codim]++;
611 int a1 = (3*k + 12)*k + 11;
613 unsigned int dof_count[16] = {0};
615 for (i[3] = 0; i[3] <= k; ++i[3])
616 for (i[2] = 0; i[2] <= k - i[3]; ++i[2])
617 for (i[1] = 0; i[1] <= k - i[2] - i[3]; ++i[1])
619 i[0] = k - i[1] - i[2] - i[3];
621 unsigned int entity = 0;
622 unsigned int codim = 0;
623 for (
unsigned int m = 0; m < 4; ++m)
625 j[m] = i[vertexMap[m]];
626 entity += !!j[m] << m;
629 int local_index = j[3]*(a1 + (a2 + j[3])*j[3])/6
630 + j[2]*(2*(k - j[3]) + 3 - j[2])/2 + j[1];
631 localKeys_[local_index] = LocalKey(subindex[entity], codim, dof_count[entity]++);
640 template<
class LocalBasis>
641 class LagrangeSimplexLocalInterpolation
643 static const int kdiv = (LocalBasis::order() == 0 ? 1 : LocalBasis::order());
653 template<
typename F,
typename C>
654 void interpolate (
const F& f, std::vector<C>& out)
const
656 constexpr auto dim = LocalBasis::Traits::dimDomain;
657 constexpr auto k = LocalBasis::order();
658 using D =
typename LocalBasis::Traits::DomainFieldType;
660 typename LocalBasis::Traits::DomainType x;
662 out.resize(LocalBasis::size());
667 auto center = ReferenceElements<D,dim>::simplex().position(0,0);
676 std::fill(x.begin(), x.end(), 0);
680 for (
int i=0; i<dim; i++)
682 for (
int j=0; j<dim; j++)
692 for (
unsigned int i=0; i<k+1; i++)
703 for (
unsigned int j=0; j<=k; j++)
704 for (
unsigned int i=0; i<=k-j; i++)
706 x = { ((D)i)/k, ((D)j)/k };
714 DUNE_THROW(NotImplemented,
"LagrangeSimplexLocalInterpolation only implemented for dim<=3!");
717 for (
int i2 = 0; i2 <= (int)k; i2++)
718 for (
int i1 = 0; i1 <= (int)k-i2; i1++)
719 for (
int i0 = 0; i0 <= (int)k-i1-i2; i0++)
721 x[0] = ((D)i0)/((D)kdiv);
722 x[1] = ((D)i1)/((D)kdiv);
723 x[2] = ((D)i2)/((D)kdiv);
788 template<
class D,
class R,
int d,
int k>
795 Impl::LagrangeSimplexLocalCoefficients<d,k>,
796 Impl::LagrangeSimplexLocalInterpolation<Impl::LagrangeSimplexLocalBasis<D,R,d,k> > >;
803 template<
typename VertexMap>
805 : coefficients_(vertexmap)
819 return coefficients_;
826 return interpolation_;
830 static constexpr std::size_t
size ()
832 return Impl::LagrangeSimplexLocalBasis<D,R,d,k>::size();
837 static constexpr GeometryType
type ()
839 return GeometryTypes::simplex(d);
843 Impl::LagrangeSimplexLocalBasis<D,R,d,k> basis_;
845 Impl::LagrangeSimplexLocalInterpolation<Impl::LagrangeSimplexLocalBasis<D,R,d,k> > interpolation_;
Definition bdfmcube.hh:18
const EdgeS0_5FiniteElement< Geometry, RF >::Traits::Coefficients & EdgeS0_5FiniteElement< Geometry, RF >::coefficients_
Definition edges0.5.hh:84
Dune::LocalBasisTraits< D, dim, FieldVector< D, dim >, R, 1, FieldVector< R, 1 >, FieldMatrix< R, 1, dim > >::DomainType D DomainType
Definition common/localbasis.hh:43
traits helper struct
Definition localfiniteelementtraits.hh:13
Impl::LagrangeSimplexLocalBasis< D, R, d, k > LocalBasisType
Definition localfiniteelementtraits.hh:16
Impl::LagrangeSimplexLocalCoefficients< d, k > LocalCoefficientsType
Definition localfiniteelementtraits.hh:20
Impl::LagrangeSimplexLocalInterpolation< Impl::LagrangeSimplexLocalBasis< D, R, d, k > > LocalInterpolationType
Definition localfiniteelementtraits.hh:24
const Traits::LocalInterpolationType & localInterpolation() const
Returns object that evaluates degrees of freedom.
Definition lagrangesimplex.hh:824
const Traits::LocalBasisType & localBasis() const
Returns the local basis, i.e., the set of shape functions.
Definition lagrangesimplex.hh:810
LagrangeSimplexLocalFiniteElement()
Definition lagrangesimplex.hh:799
LocalFiniteElementTraits< Impl::LagrangeSimplexLocalBasis< D, R, d, k >, Impl::LagrangeSimplexLocalCoefficients< d, k >, Impl::LagrangeSimplexLocalInterpolation< Impl::LagrangeSimplexLocalBasis< D, R, d, k > > > Traits
Export number types, dimensions, etc.
Definition lagrangesimplex.hh:794
static constexpr std::size_t size()
The number of shape functions.
Definition lagrangesimplex.hh:830
LagrangeSimplexLocalFiniteElement(const VertexMap &vertexmap)
Constructs a finite element given a vertex reordering.
Definition lagrangesimplex.hh:804
const Traits::LocalCoefficientsType & localCoefficients() const
Returns the assignment of the degrees of freedom to the element subentities.
Definition lagrangesimplex.hh:817
static constexpr GeometryType type()
The reference element that the local finite element is defined on.
Definition lagrangesimplex.hh:837