DUNE-FEM (unstable)

Dune::Fem::LocalFunctionAdapter< LocalFunctionImpl > Class Template Reference

LocalFunctionAdapter wrapped a class with a local evaluate method into a grid function. More...

#include <dune/fem/function/common/localfunctionadapter.hh>

Public Types

typedef LocalFunctionImpl LocalFunctionImplType
 Evaluate class.
 
typedef BaseType::FunctionType FunctionType
 type of function
 
typedef LocalFunctionAdapterTraits< LocalFunctionImplTypeTraits
 traits class
 
typedef Traits::GridPartType GridPartType
 type of grid part
 
typedef Traits::DiscreteFunctionSpaceType DiscreteFunctionSpaceType
 type of discrete function space
 
typedef Traits::GridType GridType
 type of grid
 
typedef Traits::DomainFieldType DomainFieldType
 domain type
 
typedef Traits::RangeFieldType RangeFieldType
 range type
 
typedef Traits::DomainType DomainType
 domain type
 
typedef Traits::RangeType RangeType
 range type
 
typedef Traits::JacobianRangeType JacobianRangeType
 jacobian type
 
typedef Traits::EntityType EntityType
 type of codim 0 entity
 
typedef Traits::LocalFunctionType LocalFunctionType
 type of local function to export
 
typedef LocalFunctionImpl::FunctionSpaceType FunctionSpaceType
 type of function space this function belongs to
 
typedef FunctionSpaceType ::HessianRangeType HessianRangeType
 hessian type
 
typedef Mapping< DomainFieldType, RangeFieldType, DomainType, RangeTypeMappingType
 type of mapping base class
 

Public Member Functions

 LocalFunctionAdapter (const std::string &name, LocalFunctionImplType &localFunctionImpl, const GridPartType &gridPart, unsigned int order=DiscreteFunctionSpaceType::polynomialOrder)
 
 LocalFunctionAdapter (const std::string &name, const LocalFunctionImplType &localFunctionImpl, const GridPartType &gridPart, unsigned int order=DiscreteFunctionSpaceType::polynomialOrder)
 constructor taking a const reference instance of the local function class
 
 LocalFunctionAdapter (const std::string &name, LocalFunctionImplType &&localFunctionImpl, const GridPartType &gridPart, unsigned int order=DiscreteFunctionSpaceType::polynomialOrder)
 constructor taking a r-value reference instance of the local function class
 
template<class ... Args>
 LocalFunctionAdapter (const std::string &name, const GridPartType &gridPart, unsigned int order, Args &... args)
 
template<class ... Args>
 LocalFunctionAdapter (const std::string &name, const GridPartType &gridPart, const std::tuple< Args &... > &args, unsigned int order=DiscreteFunctionSpaceType::polynomialOrder)
 
unsigned int order () const
 return the order of the space
 
bool continuous () const
 return true, probably
 
const LocalFuncStorageType & localFunctionImpl () const
 return local function implementation
 
LocalFuncStorageType & localFunctionImpl ()
 return local function implementation
 
void evaluate (const DomainType &global, RangeType &result) const
 evaluate function on local coordinate local
 
LocalFunctionType localFunction (const EntityType &entity)
 obtain a local function for an entity (read-write) More...
 
const LocalFunctionType localFunction (const EntityType &entity) const
 obtain a local function for an entity (read-write) More...
 
const std::string & name () const
 obtain the name of the discrete function More...
 
template<class DFType >
DiscreteFunctionTypeoperator+= (const DFType &g)
 
template<class DFType >
DiscreteFunctionTypeoperator-= (const DFType &g)
 substract all degrees of freedom from given discrete function using the dof iterators More...
 
DiscreteFunctionTypeoperator*= (const RangeFieldType &scalar)
 multiply all DoFs with a scalar factor More...
 
DiscreteFunctionTypeoperator/= (const RangeFieldType &scalar)
 devide all DoFs by a scalar factor More...
 
template<class ArgumentType >
void initialize (const ArgumentType &arg, double time)
 initialize local function with argument and time
 
void registerLocalFunction (LocalFunctionType *lf) const
 add LocalFunction to list of local functions
 
void deleteLocalFunction (LocalFunctionType *lf) const
 remove LocalFunction to list of local functions
 
virtual void operator() (const DomainType &arg, RangeType &dest) const
 application operator call evaluate More...
 
void evaluate (const DomainType &x, RangeType &value) const
 evaluate the function More...
 
void jacobian (const DomainType &x, JacobianRangeType &jacobian) const
 evaluate the Jacobian of the function More...
 
void hessian (const DomainType &x, HessianRangeType &hessian) const
 evaluate the hessian of the function More...
 

Related Functions

(Note that these are not member functions.)

template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator+ (const Mapping< DFieldType, RFieldType, DType, RType > &a, const Mapping< DFieldType, RFieldType, DType, RType > &b)
 add two mappings More...
 
template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator- (const Mapping< DFieldType, RFieldType, DType, RType > &a, const Mapping< DFieldType, RFieldType, DType, RType > &b)
 substract two mappings More...
 
template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator* (const Mapping< DFieldType, RFieldType, DType, RType > &mapping, const RFieldType &factor)
 scale mapping with factor More...
 
template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator* (const RFieldType &factor, const Mapping< DFieldType, RFieldType, DType, RType > &mapping)
 scale mapping with factor More...
 
template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator/ (const Mapping< DFieldType, RFieldType, DType, RType > &mapping, const RFieldType &factor)
 operator / for mappings More...
 
template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator/ (const RFieldType &factor, const Mapping< DFieldType, RFieldType, DType, RType > &mapping)
 operator / for mappings More...
 

Detailed Description

template<class LocalFunctionImpl>
class Dune::Fem::LocalFunctionAdapter< LocalFunctionImpl >

LocalFunctionAdapter wrapped a class with a local evaluate method into a grid function.

The class takes one template argument LocalFunctionImpl which holds the method evaluate(...) to evaluate the local function

template<class PointType> LocalFunctionImpl::evaluate(const PointType& x,RangeType& val)

and a method init(...)

LocalFunctionImpl::init(const EntityType& entity)

to set the entity.

It is important to know that the point x it is not necessary of type DomainType. More precisely, if the evaluate(...) is used with a caching quadrature point the type is different. Indeed floating point coordinates are not very well suited to address the cache therefore quadrature[i] return a QuadraturePointWrapper (which simply stores a reference to the quadrature and the index i).

In order to be sure that the point x is of type DomainType, you can use the function coordinate(x) which can be also called with a DomainType.

Therefore, the local implementation should be something like

template<class PointType> LocalFunctionImpl::evaluate(const PointType& x,RangeType& val) { const DomainType xDomain(coordiante(x)); // do stuff with xDomain }

to avoid type conflicts.

Required type in LocalFunctionImpl are:

FunctionSpaceType GridPartType EntityType DomainType RangeType

An instance of the LocalFunctionImpl class is passed to the constructor.

In order to adapt a lambda or a plain C++ function, you can directly use the LocalAnalyticalFunctionBinder which provides all the necessary types and methods.

Member Function Documentation

◆ evaluate()

void Dune::Fem::Function< LocalFunctionImpl::FunctionSpaceType , LocalFunctionAdapter< LocalFunctionImpl > >::evaluate ( const DomainType x,
RangeType value 
) const
inlineinherited

evaluate the function

Parameters
[in]xevaluation point
[out]valuevalue of the function in x

◆ hessian()

void Dune::Fem::Function< LocalFunctionImpl::FunctionSpaceType , LocalFunctionAdapter< LocalFunctionImpl > >::hessian ( const DomainType x,
HessianRangeType hessian 
) const
inlineinherited

evaluate the hessian of the function

Parameters
[in]xevaluation point
[out]hessianvalue of the hessian in x

◆ jacobian()

void Dune::Fem::Function< LocalFunctionImpl::FunctionSpaceType , LocalFunctionAdapter< LocalFunctionImpl > >::jacobian ( const DomainType x,
JacobianRangeType jacobian 
) const
inlineinherited

evaluate the Jacobian of the function

Parameters
[in]xevaluation point
[out]jacobianvalue of the Jacobian in x

◆ operator()()

virtual void Dune::Fem::Function< LocalFunctionImpl::FunctionSpaceType , LocalFunctionAdapter< LocalFunctionImpl > >::operator() ( const DomainType arg,
RangeType dest 
) const
inlinevirtualinherited

application operator call evaluate

Parameters
[in]argargument
[out]destdestination, i.e. f(arg)

Friends And Related Function Documentation

◆ operator*() [1/2]

template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator* ( const Mapping< DFieldType, RFieldType, DType, RType > &  mapping,
const RFieldType &  factor 
)
related

scale mapping with factor

Parameters
[in]mappingMapping which is scaled
[in]factorfactor with which mapping is scaled
Returns
new object mapping

References Dune::Fem::MappingOperators::multiplyMapping().

◆ operator*() [2/2]

template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator* ( const RFieldType &  factor,
const Mapping< DFieldType, RFieldType, DType, RType > &  mapping 
)
related

scale mapping with factor

Parameters
[in]factorfactor with which mapping is scaled
[in]mappingMapping which is scaled
Returns
new object mapping

References Dune::Fem::MappingOperators::multiplyMapping().

◆ operator+()

template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator+ ( const Mapping< DFieldType, RFieldType, DType, RType > &  a,
const Mapping< DFieldType, RFieldType, DType, RType > &  b 
)
related

add two mappings

Parameters
[in]amapping 1
[in]bmapping 2
Returns
new object mapping

References Dune::Fem::MappingOperators::addMappings().

◆ operator-()

template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator- ( const Mapping< DFieldType, RFieldType, DType, RType > &  a,
const Mapping< DFieldType, RFieldType, DType, RType > &  b 
)
related

substract two mappings

Parameters
[in]amapping 1
[in]bmapping 2
Returns
new object mapping

References Dune::Fem::MappingOperators::substractMappings().

◆ operator/() [1/2]

template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator/ ( const Mapping< DFieldType, RFieldType, DType, RType > &  mapping,
const RFieldType &  factor 
)
related

operator / for mappings

Parameters
[in]mappingmapping which is divided
[in]factorf factor by which result of mapping is divided
Returns
new object mapping

References Dune::Fem::MappingOperators::divideMapping().

◆ operator/() [2/2]

template<class DFieldType , class RFieldType , class DType , class RType >
static Mapping< DFieldType, RFieldType, DType, RType > operator/ ( const RFieldType &  factor,
const Mapping< DFieldType, RFieldType, DType, RType > &  mapping 
)
related

operator / for mappings

Parameters
[in]factorby which result of mapping is divided
[in]mappingwhich is divided
Returns
new object mapping

References Dune::Fem::MappingOperators::divideMapping().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 4, 22:30, 2024)