6#ifndef DUNE_PDELAB_FINITEELEMENT_QKDGLEGENDRE_HH
7#define DUNE_PDELAB_FINITEELEMENT_QKDGLEGENDRE_HH
18#include <dune/localfunctions/common/localbasis.hh>
19#include <dune/localfunctions/common/localfiniteelementtraits.hh>
20#include <dune/localfunctions/common/localinterpolation.hh>
21#include <dune/localfunctions/common/localkey.hh>
22#include <dune/localfunctions/common/localtoglobaladaptors.hh>
23#include <dune/localfunctions/common/localinterpolation.hh>
28 namespace LegendreStuff
33 template<
int k,
int n>
37 value=(k+1)*LegendreSize<k,n-1>::value
42 struct LegendreSize<0,1>
50 struct LegendreSize<k,1>
58 struct LegendreSize<0,n>
65 template<
int k,
int d>
69 for (
int j=0; j<d; j++)
78 template<
class D,
class R,
int k>
83 void p (D x, std::vector<R>& value)
const
88 for (
int n=2; n<=k; n++)
89 value[n] = ((2*n-1)*(2*x-1)*value[n-1]-(n-1)*value[n-2])/n;
93 R
p (
int i, D x)
const
95 std::vector<R> v(k+1);
101 void dp (D x, std::vector<R>&
derivative)
const
109 for (
int n=2; n<=k; n++)
111 value[n] = ((2*n-1)*(2*x-1)*value[n-1]-(n-1)*value[n-2])/n;
117 void pdp (D x, std::vector<R>& value, std::vector<R>&
derivative)
const
125 for (
int n=2; n<=k; n++)
127 value[n] = ((2*n-1)*(2*x-1)*value[n-1]-(n-1)*value[n-2])/n;
133 R dp (
int i, D x)
const
135 std::vector<R> v(k+1);
141 template<
class D,
class R>
142 class LegendrePolynomials1d<D,R,0>
146 void p (D x, std::vector<R>& value)
const
153 R
p (
int i, D x)
const
159 void dp (D x, std::vector<R>&
derivative)
const
166 R dp (
int i, D x)
const
172 void pdp (D x, std::vector<R>& value, std::vector<R>&
derivative)
const
181 template<
class D,
class R>
182 class LegendrePolynomials1d<D,R,1>
186 void p (D x, std::vector<R>& value)
const
194 R
p (
int i, D x)
const
196 return (1-i) + i*(2*x-1);
200 void dp (D x, std::vector<R>&
derivative)
const
208 R dp (
int i, D x)
const
210 return (1-i)*0 + i*(2);
214 void pdp (D x, std::vector<R>& value, std::vector<R>&
derivative)
const
237 template<
class D,
class R,
int k,
int d>
240 enum { n = LegendreSize<k,d>::value };
242 mutable std::vector<std::vector<R> > v;
243 mutable std::vector<std::vector<R> > a;
246 typedef LocalBasisTraits<D,d,Dune::FieldVector<D,d>,R,1,
Dune::FieldVector<R,1>,
Dune::FieldMatrix<R,1,d> >
Traits;
259 std::vector<typename Traits::RangeType>& value)
const
265 for (
size_t j=0; j<d; j++) poly.p(x[j],v[j]);
268 for (
size_t i=0; i<n; i++)
277 for (
int j=0; j<d; j++) value[i] *= v[j][alpha[j]];
284 std::vector<typename Traits::JacobianType>& value)
const
287 value.resize(
size());
290 for (
size_t j=0; j<d; j++) poly.pdp(x[j],v[j],a[j]);
293 for (
size_t i=0; i<n; i++)
299 for (
int j=0; j<d; j++)
302 value[i][0][j] = a[j][alpha[j]];
305 for (
int l=0; l<d; l++)
307 value[i][0][j] *= v[l][alpha[l]];
315 std::vector<typename Traits::RangeType>& out)
const {
317 if (totalOrder == 0) {
333 template<
class D,
class R,
int k,
int d>
336 enum { n = LegendreSize<k,d>::value };
347 template<
typename F,
typename C>
353 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
359 std::vector<R> diagonal(n);
360 for (
int i=0; i<n; i++) { out[i] = 0.0; diagonal[i] = 0.0; }
364 it=rule.begin(); it!=rule.end(); ++it)
369 for (
int i=0; i<d; i++) x[i] = it->position()[i];
373 std::vector<RangeType> phi(n);
377 for (
int i=0; i<n; i++) {
378 out[i] += y*phi[i]*it->weight();
379 diagonal[i] += phi[i]*phi[i]*it->weight();
382 for (
int i=0; i<n; i++) out[i] /= diagonal[i];
392 template<
int k,
int d>
395 enum { n = LegendreSize<k,d>::value };
401 for (std::size_t i=0; i<n; i++)
418 std::vector<LocalKey> li;
425 template<
class D,
class R,
int k,
int d>
434 enum { n = LegendreStuff::LegendreSize<k,d>::value };
443 : gt(GeometryTypes::
cube(d))
464 return interpolation;
488 LocalCoefficients coefficients;
489 LocalInterpolation interpolation;
499 template<
class Geometry,
class RF,
int k>
502 QkDGLegendreLocalFiniteElement<
503 typename Geometry::ctype, RF, k, Geometry::mydimension
513 static const LFE lfe;
520 template<
class Geometry,
class RF,
int k>
521 const typename QkDGLegendreFiniteElementFactory<Geometry, RF, k>::LFE
522 QkDGLegendreFiniteElementFactory<Geometry, RF, k>::lfe;
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:92
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:100
static constexpr int mydimension
geometry dimension
Definition: geometry.hh:94
Lagrange shape functions of order k on the reference cube.
Definition: qkdglegendre.hh:239
unsigned int size() const
number of shape functions
Definition: qkdglegendre.hh:252
void evaluateFunction(const typename Traits::DomainType &x, std::vector< typename Traits::RangeType > &value) const
Evaluate all shape functions.
Definition: qkdglegendre.hh:258
void evaluateJacobian(const typename Traits::DomainType &x, std::vector< typename Traits::JacobianType > &value) const
Evaluate Jacobian of all shape functions.
Definition: qkdglegendre.hh:283
void partial(const std::array< unsigned int, Traits::dimDomain > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivative of all shape functions.
Definition: qkdglegendre.hh:313
unsigned int order() const
Polynomial order of the shape functions.
Definition: qkdglegendre.hh:325
Layout map for Q1 elements.
Definition: qkdglegendre.hh:394
DGLegendreLocalCoefficients()
Standard constructor.
Definition: qkdglegendre.hh:399
std::size_t size() const
number of coefficients
Definition: qkdglegendre.hh:406
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: qkdglegendre.hh:412
determine degrees of freedom
Definition: qkdglegendre.hh:335
void interpolate(const F &ff, std::vector< C > &out) const
Local interpolation of a function.
Definition: qkdglegendre.hh:348
The 1d Legendre Polynomials (k=0,1 are specialized below)
Definition: qkdglegendre.hh:80
void p(D x, std::vector< R > &value) const
evaluate all polynomials at point x
Definition: qkdglegendre.hh:83
Describe position of one degree of freedom.
Definition: localkey.hh:24
Default exception for dummy implementations.
Definition: exceptions.hh:355
Factory for global-valued DGLegendre elements.
Definition: qkdglegendre.hh:507
QkDGLegendreFiniteElementFactory()
default constructor
Definition: qkdglegendre.hh:517
Definition: qkdglegendre.hh:427
const Traits::LocalCoefficientsType & localCoefficients() const
Definition: qkdglegendre.hh:455
LocalFiniteElementTraits< LocalBasis, LocalCoefficients, LocalInterpolation > Traits
Definition: qkdglegendre.hh:438
const Traits::LocalInterpolationType & localInterpolation() const
Definition: qkdglegendre.hh:462
GeometryType type() const
Definition: qkdglegendre.hh:476
std::size_t size() const
Definition: qkdglegendre.hh:469
const Traits::LocalBasisType & localBasis() const
Definition: qkdglegendre.hh:448
QkDGLegendreLocalFiniteElement()
Definition: qkdglegendre.hh:442
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:214
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:326
Factory for ScalarLocalToGlobalFiniteElementAdaptor objects.
Definition: localtoglobaladaptors.hh:244
Definition of the DUNE_NO_DEPRECATED_* macros.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E,...)
Definition: exceptions.hh:312
TrigonometricFunction< K, -cosFactor, sinFactor > derivative(const TrigonometricFunction< K, sinFactor, cosFactor > &f)
Obtain derivative of TrigonometricFunction function.
Definition: trigonometricfunction.hh:43
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:462
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:280
Dune namespace.
Definition: alignedallocator.hh:13
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:35
D DomainType
domain type
Definition: localbasis.hh:43
R RangeType
range type
Definition: localbasis.hh:52
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
A unique label for each type of element that can occur in a grid.