1 #ifndef DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH 2 #define DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH 7 #include <dune/grid/common/capabilities.hh> 22 template<
class GridPart,
class BasisFunctionSet,
23 class Quadrature = CachingQuadrature< GridPart, BasisFunctionSet::EntityType::codimension > >
31 template<
class Gr
idPart,
class BasisFunctionSet,
class Quadrature >
33 :
public LocalL2Projection< BasisFunctionSet, DiscontinuousGalerkinLocalL2Projection< GridPart, BasisFunctionSet, Quadrature > >
43 typedef GridPart GridPartType;
44 typedef typename GridPartType::GridType GridType;
45 static const bool cartesian = Dune::Capabilities::isCartesian< GridType >::v;
47 typedef typename std::conditional< cartesian,
50 >::type LocalRieszProjectionType;
60 : impl_( LocalRieszProjectionType( basisFunctionSet ) )
64 : impl_( LocalRieszProjectionType(
std::forward< BasisFunctionSetType >(
basisFunctionSet ) ) )
76 : impl_(
std::
move( other.impl_ ) )
79 ThisType &
operator= (
const ThisType & ) =
default;
100 template<
class LocalFunction,
class LocalDofVector >
103 impl_( localFunction, localDofVector );
109 Implementation impl_;
116 #endif // #ifndef DUNE_FEM_SPACE_DISCONTINUOUSGALERKIN_INTERPOLATION_HH
BaseType::BasisFunctionSetType BasisFunctionSetType
basis function set type
Definition: discontinuousgalerkin/interpolation.hh:40
DiscontinuousGalerkinLocalL2Projection(BasisFunctionSetType &&basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:63
void apply(const LocalFunction &localFunction, LocalDofVector &localDofVector) const
please doc me
Definition: discontinuousgalerkin/interpolation.hh:101
interface for local functions
Definition: localfunction.hh:41
Definition: coordinate.hh:4
DiscontinuousGalerkinLocalL2Projection(const BasisFunctionSetType &basisFunctionSet)
Definition: discontinuousgalerkin/interpolation.hh:59
BasisFunctionSet basisFunctionSet() const
return basis function set
Definition: discontinuousgalerkin/interpolation.hh:94
Definition: operator/projection/local/riesz/orthonormal.hh:23
Definition: discontinuousgalerkin/interpolation.hh:24
void move(ArrayInterface< T > &array, const unsigned int oldOffset, const unsigned int newOffset, const unsigned int length)
Definition: array_inline.hh:38
please doc me
Definition: local/l2projection.hh:27
ThisType & operator=(const ThisType &)=default
BasisFunctionSetType basisFunctionSet() const
return basis function set
Definition: local/l2projection.hh:175
DiscontinuousGalerkinLocalL2Projection(ThisType &&other)
Definition: discontinuousgalerkin/interpolation.hh:75
Definition: basisfunctionset/basisfunctionset.hh:31