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 )
110 template<
class Po
int >
111 void evaluate (
const Point &p, RangeType &ret )
const 118 template<
class Po
int >
119 void jacobian (
const Point &p, JacobianRangeType &jac )
const 121 HostJacobianRangeType hJac;
126 template<
class Po
int >
127 void hessian (
const Point &p, HessianRangeType &hes )
const 129 HostHessianRangeType hHes;
130 this->
get().
hessian( p, hHes );
134 template<
class Quadrature,
class Vector >
147 template<
class QuadratureType,
class VectorType >
150 const unsigned int nop = quadrature.nop();
151 for(
unsigned int qp = 0; qp < nop; ++qp )
152 evaluate( quadrature[ qp ], values[ qp ] );
155 template<
class QuadratureType,
class VectorType >
156 void evaluateQuadratureImp (
const QuadratureType &quadrature, VectorType &values,
const JacobianRangeType & )
const 158 const unsigned int nop = quadrature.nop();
159 for(
unsigned int qp = 0; qp < nop; ++qp )
160 jacobian( quadrature[ qp ], values[ qp ] );
171 template<
class HostLocalFunction,
class Converter >
176 return LocalFunctionConverterType(
std::move( hostLocalFunction ), converter );
179 template<
class HostLocalFunction,
class Converter >
181 localFunctionConverter ( std::reference_wrapper< HostLocalFunction > hostLocalFunction,
const Converter &converter = Converter() )
184 return LocalFunctionConverterType( hostLocalFunction.get(), converter );
191 #endif // #ifndef DUNE_FEM_FUNCTION_LOCALFUNCTION_CONVERTER_HH
LocalFunctionConverter(const HostLocalFunction &hostLocalFunction, const Converter &converter=Converter())
Definition: converter.hh:102
FunctionSpaceType::RangeFieldType RangeFieldType
Definition: converter.hh:98
void jacobian(const Point &p, JacobianRangeType &jac) const
Definition: converter.hh:119
FunctionSpaceType::RangeType RangeType
Definition: converter.hh:94
LocalFunctionConverter(HostLocalFunction &&hostLocalFunction, const Converter &converter=Converter())
Definition: converter.hh:106
static const int dimDomain
Definition: converter.hh:100
HostLocalFunction::EntityType EntityType
Definition: converter.hh:90
void evaluate(const Point &p, RangeType &ret) const
Definition: converter.hh:111
Converter converter_
Definition: converter.hh:163
void hessian(const Point &p, HessianRangeType &hes) const
Definition: converter.hh:127
ToNewDimRangeFunctionSpace< typename HostLocalFunction::FunctionSpaceType, dimRange >::Type FunctionSpaceType
Definition: converter.hh:87
static const int dimRange
Definition: converter.hh:84
int order() const
Definition: converter.hh:140
Definition: coordinate.hh:4
LocalFunctionConverter< HostLocalFunction, Converter, __InstationaryFunction::HoldCopy > localFunctionConverter(HostLocalFunction hostLocalFunction, const Converter &converter=Converter())
Definition: converter.hh:173
void evaluateQuadrature(const Quadrature &quad, Vector &vector) const
Definition: converter.hh:135
FunctionSpaceType::JacobianRangeType JacobianRangeType
Definition: converter.hh:95
FunctionSpaceType::DomainFieldType DomainFieldType
Definition: converter.hh:97
void evaluateQuadratureImp(const QuadratureType &quadrature, VectorType &values, const JacobianRangeType &) const
Definition: converter.hh:156
implementation of a Dune::Fem::LocalFunction on a FunctionSpace V restircted/prolongated from an othe...
Definition: converter.hh:72
void move(ArrayInterface< T > &array, const unsigned int oldOffset, const unsigned int newOffset, const unsigned int length)
Definition: array_inline.hh:38
void evaluateQuadratureImp(const QuadratureType &quadrature, VectorType &values, const RangeType &) const
Definition: converter.hh:148
FunctionSpaceType::DomainType DomainType
Definition: converter.hh:93
const EntityType & entity() const
Definition: converter.hh:142
void init(const EntityType &entity)
Definition: converter.hh:144
FunctionSpaceType::HessianRangeType HessianRangeType
Definition: converter.hh:96
convert functions space to space with new dim range
Definition: functionspace.hh:246
actual interface class for quadratures
Definition: quadrature.hh:320