DUNE-ACFEM (2.5.1)

If a GridFunctionExpression only consists of wrapped components, constants and parameters, then it is much more efficient to wrap the entire expression once more in order to make sure that the transformation from the local to the global coordinate system only has to be performed once for the entire expression instead of computing the global coordinates again and again for each component of the expression. More...

Classes

class  Dune::ACFem::UnaryGridFunctionExpression< UnOp, FunctionType >
 UnaryGridFunctionExpression implements point-wise unary operations for GridFunction-types. More...
 
class  Dune::ACFem::BinaryGridFunctionExpression< BinOp, LeftFunction, RightFunction >
 BinaryGridFunctionExpression implements point-wise vector-space operations for GridFunction-types. More...
 

Functions

template<class Function , class GridPart >
static auto Dune::ACFem::operator- (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(-asExprArg(f_)), GridPart >
 -Wrapped(f) = Wrapped(-f)
 
template<class Function , class GridPart >
static auto Dune::ACFem::operator+ (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< Function, GridPart >
 +Wrapped(f) = Wrapped(f)
 
template<class Left , class Right , class GridPart >
static auto Dune::ACFem::operator+ (const GridFunctionWrapper< Left, GridPart > &f_, const GridFunctionWrapper< Right, GridPart > &g_) -> GridFunctionWrapper< decltype(asExprArg(f_)+asExprArg(g_)), GridPart >
 Wrapped(f) + Wrapped(g) = Wrapped(f + g)
 
template<class Left , class Right , class GridPart >
static auto Dune::ACFem::operator- (const GridFunctionWrapper< Left, GridPart > &f_, const GridFunctionWrapper< Right, GridPart > &g_) -> GridFunctionWrapper< decltype(asExprArg(f_) - asExprArg(g_)), GridPart >
 Wrapped(f) - Wrapped(g) = Wrapped(f - g)
 
template<class Function , class GridPart >
static auto Dune::ACFem::operator* (const typename Function::RangeFieldType &s_, const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(s_ *asExprArg(f_)), GridPart >
 s * Wrapped(f) = Wrapped(s * f)
 
template<class Function , class GridPart >
static auto Dune::ACFem::operator* (const GridFunctionWrapper< Function, GridPart > &f_, const typename Function::RangeFieldType &s_) -> decltype(s_ *f_)
 Wrapped(f) * s = Wrapped(s * f)
 
template<class Parameter , class Function , class GridPart >
static auto Dune::ACFem::operator* (const ParameterInterface< Parameter > &p_, const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(asImp(p_) *asExprArg(f_)), GridPart >
 Param * Wrapped(f) = Wrapped(Param * f)
 
template<class Parameter , class Function , class GridPart >
static auto Dune::ACFem::operator* (const GridFunctionWrapper< Function, GridPart > &f_, const ParameterInterface< Parameter > &p_) -> decltype(asImp(p_) *f_)
 Wrapped(f) * Param = Wrapped(Param * f)
 
template<class Field , class Parameter , class Function , class GridPart >
static auto Dune::ACFem::operator* (const BinaryParameterExpression< SMultiplyOperation, Field, Parameter > &p_, const GridFunctionWrapper< Function, GridPart > &f_) -> decltype(p_.left() *(p_.right() *f_))
 Move scalars out of BinaryParameterExpression.
 
template<class Left , class Right , class GridPart >
static auto Dune::ACFem::operator* (const GridFunctionWrapper< Left, GridPart > &f_, const GridFunctionWrapper< Right, GridPart > &g_) -> GridFunctionWrapper< decltype(asExprArg(f_) *asExprArg(g_)), GridPart >
 Wrapped(f) * Wrapped(g) = Wrapped(f * g)
 
template<class Left , class Right , class GridPart >
static auto Dune::ACFem::operator/ (const GridFunctionWrapper< Left, GridPart > &f_, const GridFunctionWrapper< Right, GridPart > &g_) -> GridFunctionWrapper< decltype(asExprArg(f_)/asExprArg(g_)), GridPart >
 Wrapped(f) / Wrapped(scalar g) = Wrapped(f * 1/g)
 
template<class Function , class GridPart >
static auto Dune::ACFem::operator/ (const typename Function::RangeFieldType &s_, const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(s_/asExprArg(f_)), GridPart >
 s / Wrapped(g) = Wrapped(s / g)
 
template<class Parameter , class Function , class GridPart >
static auto Dune::ACFem::operator/ (const ParameterInterface< Parameter > &p_, const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(asImp(p_)/asExprArg(f_)), GridPart >
 Param / Wrapped(f) = Wrapped(Param / f)
 
template<class Function , class GridPart >
static auto Dune::ACFem::exp (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(exp(asExprArg(f_))), GridPart >
 exp(Wrapped(f)) = Wrapped(exp(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::log (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(log(asExprArg(f_))), GridPart >
 log(Wrapped(f)) = Wrapped(log(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::sqrt (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(sqrt(asExprArg(f_))), GridPart >
 sqrt(Wrapped(f)) = Wrapped(sqrt(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::sqr (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(sqr(asExprArg(f_))), GridPart >
 sqr(Wrapped(f)) = Wrapped(sqr(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::sin (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(sin(asExprArg(f_))), GridPart >
 sin(Wrapped(f)) = Wrapped(sin(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::cos (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(cos(asExprArg(f_))), GridPart >
 cos(Wrapped(f)) = Wrapped(cos(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::tan (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(tan(asExprArg(f_))), GridPart >
 tan(Wrapped(f)) = Wrapped(tan(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::atan (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(atan(asExprArg(f_))), GridPart >
 atan(Wrapped(f)) = Wrapped(atan(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::asin (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(asin(asExprArg(f_))), GridPart >
 asin(Wrapped(f)) = Wrapped(asin(f))
 
template<class Function , class GridPart >
static auto Dune::ACFem::acos (const GridFunctionWrapper< Function, GridPart > &f_) -> GridFunctionWrapper< decltype(acos(asExprArg(f_))), GridPart >
 acos(Wrapped(f)) = Wrapped(acos(f))
 

Detailed Description

If a GridFunctionExpression only consists of wrapped components, constants and parameters, then it is much more efficient to wrap the entire expression once more in order to make sure that the transformation from the local to the global coordinate system only has to be performed once for the entire expression instead of computing the global coordinates again and again for each component of the expression.

This works similar to adding BoundaryIndicators in the BoundaryFunctionExpressions. We do not remove the inner wrappers in order to make us of the GridFunctionExpression templates.

Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 24, 22:30, 2024)