1#ifndef DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH
2#define DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH
7#include <dune/grid/common/capabilities.hh>
24 template<
class GridPart,
class BasisFunctionSet,
25 class Quadrature = CachingQuadrature< GridPart, BasisFunctionSet::EntityType::codimension > >
26 class DiscontinuousGalerkinLocalL2Projection;
28 template<
class GridPart,
class BasisFunctionSet,
29 class Quadrature = CachingQuadrature< GridPart, BasisFunctionSet::EntityType::codimension > >
30 class LocalOrthonormalL2Projection;
37 template<
class Gr
idPart,
class BasisFunctionSet,
class Quadrature >
39 :
public LocalL2Projection< BasisFunctionSet, LocalOrthonormalL2Projection< GridPart, BasisFunctionSet, Quadrature > >
48 typedef typename BasisFunctionSetType :: EntityType
EntityType;
51 typedef GridPart GridPartType;
52 typedef typename GridPartType::GridType GridType;
53 typedef typename BasisFunctionSetType :: RangeType RangeType;
95 template<
class LocalFunction,
class LocalDofVector >
101 if( entity.type().isNone() )
118 template <
class QuadImpl,
class LocalFunction,
class LocalDofVector >
120 const QuadImpl& quadrature,
124 localDofVector.clear();
126 const int nop = quadrature.nop();
135 for(
auto qp : quadrature )
136 values_[ qp.index() ] *= qp.weight();
150 template<
class Gr
idPart,
class BasisFunctionSet,
class Quadrature >
152 :
public LocalL2Projection< BasisFunctionSet, DiscontinuousGalerkinLocalL2Projection< GridPart, BasisFunctionSet, Quadrature > >
162 typedef GridPart GridPartType;
163 typedef typename GridPartType::GridType GridType;
164 static const bool cartesian = Dune::Capabilities::isCartesian< GridType >::v;
166 typedef typename std::conditional< cartesian,
169 >::type LocalRieszProjectionType;
195 : impl_(
std::move( other.impl_ ) )
202 impl_ = std::move( other.impl_ );
219 template<
class LocalFunction,
class LocalDofVector >
230 Implementation impl_;
Definition: bindguard.hh:11
static GridFunctionView< GF > localFunction(const GF &gf)
Definition: gridfunctionview.hh:118
interface for local functions
Definition: localfunction.hh:77
please doc me
Definition: l2projection.hh:29
BasisFunctionSetType basisFunctionSet() const
return basis function set
Definition: l2projection.hh:179
Definition: operator/projection/local/riesz/orthonormal.hh:25
quadrature on the codim-0 reference element
Definition: elementquadrature.hh:58
actual interface class for quadratures
Definition: quadrature.hh:405
Interface class for basis function sets.
Definition: basisfunctionset/basisfunctionset.hh:31
void axpy(const Quadrature &quad, const Vector &values, DofVector &dofs) const
evaluate all basis function and multiply with given values and add to dofs
Definition: discontinuousgalerkin/interpolation.hh:153
BasisFunctionSet basisFunctionSet() const
return basis function set
Definition: discontinuousgalerkin/interpolation.hh:213
void unbind()
Definition: discontinuousgalerkin/interpolation.hh:225
BaseType::BasisFunctionSetType BasisFunctionSetType
basis function set type
Definition: discontinuousgalerkin/interpolation.hh:159
ThisType & operator=(const ThisType &)=default
DiscontinuousGalerkinLocalL2Projection(const BasisFunctionSetType &basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:178
DiscontinuousGalerkinLocalL2Projection(const ThisType &)=default
DiscontinuousGalerkinLocalL2Projection(BasisFunctionSetType &&basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:182
DiscontinuousGalerkinLocalL2Projection(ThisType &&other)
Definition: discontinuousgalerkin/interpolation.hh:194
void apply(const LocalFunction &localFunction, LocalDofVector &localDofVector) const
please doc me
Definition: discontinuousgalerkin/interpolation.hh:220
specialization of local L2 projection for orthonormal DG spaces
Definition: discontinuousgalerkin/interpolation.hh:40
BasisFunctionSetType basisFunctionSet_
Definition: discontinuousgalerkin/interpolation.hh:142
LocalOrthonormalL2Projection(ThisType &&other)=default
ThisType & operator=(const ThisType &)=default
void apply(const LocalFunction &localFunction, LocalDofVector &localDofVector) const
please doc me
Definition: discontinuousgalerkin/interpolation.hh:96
std::vector< RangeType > values_
Definition: discontinuousgalerkin/interpolation.hh:143
BasisFunctionSetType::EntityType EntityType
Definition: discontinuousgalerkin/interpolation.hh:48
LocalOrthonormalL2Projection(BasisFunctionSetType &&basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:64
void computeL2Projection(const EntityType &entity, const QuadImpl &quadrature, const LocalFunction &localFunction, LocalDofVector &localDofVector) const
Definition: discontinuousgalerkin/interpolation.hh:119
LocalOrthonormalL2Projection(const ThisType &)=default
BaseType::BasisFunctionSetType BasisFunctionSetType
basis function set type
Definition: discontinuousgalerkin/interpolation.hh:46
const BasisFunctionSet & basisFunctionSet() const
return basis function set
Definition: discontinuousgalerkin/interpolation.hh:89
LocalOrthonormalL2Projection(const BasisFunctionSetType &basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:60