1 #ifndef DUNE_FEM_FUNCTION_LOCALFUNCTION_CONVERTER_HH
2 #define DUNE_FEM_FUNCTION_LOCALFUNCTION_CONVERTER_HH
71 template<
class HostLocalFunction,
class Converter,
template<
class >
class Storage = __InstationaryFunction::HoldCopy >
73 :
private Storage< HostLocalFunction >
76 typedef Storage< HostLocalFunction > BaseType;
78 typedef typename HostLocalFunction::RangeType HostRangeType;
79 typedef typename HostLocalFunction::JacobianRangeType HostJacobianRangeType;
80 typedef typename HostLocalFunction::HessianRangeType HostHessianRangeType;
84 static const int dimRange = decltype( std::declval< Converter >() ( std::declval< HostRangeType >() ) ) ::dimension;
90 typedef typename HostLocalFunction::EntityType
EntityType;
93 typedef typename FunctionSpaceType::DomainType
DomainType;
94 typedef typename FunctionSpaceType::RangeType
RangeType;
100 static const int dimDomain = FunctionSpaceType::dimDomain;
103 : BaseType( hostLocalFunction ),
converter_( converter )
107 : BaseType( std::move( hostLocalFunction ) ),
converter_( converter )
110 template<
class Po
int >
114 this->
get().evaluate( p, hRet );
118 template<
class Po
int >
121 HostJacobianRangeType hJac;
122 this->
get().jacobian( p, hJac );
126 template<
class Po
int >
129 HostHessianRangeType hHes;
130 this->
get().hessian( p, hHes );
134 template<
class Quadrature,
class ... Vectors >
137 std::ignore = std::make_tuple(
143 const EntityType &
entity ()
const {
return this->
get().entity(); }
148 template<
class QuadratureType,
class VectorType >
151 const unsigned int nop = quadrature.nop();
152 for(
unsigned int qp = 0; qp < nop; ++qp )
153 evaluate( quadrature[ qp ], values[ qp ] );
156 template<
class QuadratureType,
class VectorType >
159 const unsigned int nop = quadrature.nop();
160 for(
unsigned int qp = 0; qp < nop; ++qp )
161 jacobian( quadrature[ qp ], values[ qp ] );
164 template<
class QuadratureType,
class VectorType >
167 const unsigned int nop = quadrature.nop();
168 for(
unsigned int qp = 0; qp < nop; ++qp )
169 hessian( quadrature[ qp ], values[ qp ] );
180 template<
class HostLocalFunction,
class Converter >
185 return LocalFunctionConverterType( std::move( hostLocalFunction ), converter );
188 template<
class HostLocalFunction,
class Converter >
189 LocalFunctionConverter< typename std::remove_const< HostLocalFunction >::type, Converter, __InstationaryFunction::HoldReference >
190 localFunctionConverter ( std::reference_wrapper< HostLocalFunction > hostLocalFunction,
const Converter &converter = Converter() )
193 return LocalFunctionConverterType( hostLocalFunction.get(), converter );
Definition: bindguard.hh:11
std::tuple_element< i, Tuple >::type & get(Dune::TypeIndexedTuple< Tuple, Types > &tuple)
Definition: typeindexedtuple.hh:122
LocalFunctionConverter< HostLocalFunction, Converter, __InstationaryFunction::HoldCopy > localFunctionConverter(HostLocalFunction hostLocalFunction, const Converter &converter=Converter())
Definition: converter.hh:182
implementation of a Dune::Fem::LocalFunction on a FunctionSpace V restircted/prolongated from an othe...
Definition: converter.hh:74
static const int dimRange
Definition: converter.hh:84
FunctionSpaceType::RangeType RangeType
Definition: converter.hh:94
FunctionSpaceType::RangeFieldType RangeFieldType
Definition: converter.hh:98
FunctionSpaceType::DomainType DomainType
Definition: converter.hh:93
void evaluateQuadrature(const Quadrature &quad, Vectors &... vector) const
Definition: converter.hh:135
FunctionSpaceType::DomainFieldType DomainFieldType
Definition: converter.hh:97
void init(const EntityType &entity)
Definition: converter.hh:145
HostLocalFunction::EntityType EntityType
Definition: converter.hh:90
int order() const
Definition: converter.hh:141
void evaluateQuadratureImp(const QuadratureType &quadrature, VectorType &values, const RangeType &) const
Definition: converter.hh:149
void evaluateQuadratureImp(const QuadratureType &quadrature, VectorType &values, const HessianRangeType &) const
Definition: converter.hh:165
void evaluateQuadratureImp(const QuadratureType &quadrature, VectorType &values, const JacobianRangeType &) const
Definition: converter.hh:157
void hessian(const Point &p, HessianRangeType &hes) const
Definition: converter.hh:127
FunctionSpaceType::JacobianRangeType JacobianRangeType
Definition: converter.hh:95
const EntityType & entity() const
Definition: converter.hh:143
Converter converter_
Definition: converter.hh:172
FunctionSpaceType::HessianRangeType HessianRangeType
Definition: converter.hh:96
LocalFunctionConverter(HostLocalFunction &&hostLocalFunction, const Converter &converter=Converter())
Definition: converter.hh:106
void jacobian(const Point &p, JacobianRangeType &jac) const
Definition: converter.hh:119
static const int dimDomain
Definition: converter.hh:100
ToNewDimRangeFunctionSpace< typename HostLocalFunction::FunctionSpaceType, dimRange >::Type FunctionSpaceType
Definition: converter.hh:87
LocalFunctionConverter(const HostLocalFunction &hostLocalFunction, const Converter &converter=Converter())
Definition: converter.hh:102
void evaluate(const Point &p, RangeType &ret) const
Definition: converter.hh:111
actual interface class for quadratures
Definition: quadrature/quadrature.hh:366
convert functions space to space with new dim range
Definition: functionspace.hh:250