KASKADE 7 development version
Modules | Functions
Working with finite element functions

Tools for representing and working with finite element functions. More...

Modules

 Function interpolation and transfer
 Tools for converting finite element functions.
 

Functions

template<class Function , class Functor >
void Kaskade::forEachQuadPoint (Function const &f, Functor &accumulator)
 Loops over the whole domain occupied by the grid and calls an accumulator for each integration point with the value of the function f. This version of forEach requires the function f to be evaluatable by an Evaluator of the associated space of f. More...
 
template<class Functor , class Space >
void Kaskade::forEachQP (Functor &g, Space const &space)
 Loops over the whole domain occupied by the grid and calls the functor g for each integration point with an evaluator of the associated space. More...
 
template<class Functor , class Function , class Space >
void Kaskade::forEachQuadPoint (Function const &f, Functor &g, Space &space)
 Loops over the whole domain occupied by the grid and calls the functor g for each integration point with the value of the function f. More...
 
template<class FEFunction >
FEFunction::ValueType Kaskade::integrateOverIntersection (FEFunction const &function, typename FEFunction::Space::GridView::Intersection const &intersection, typename FEFunction::Space::Evaluator &evaluator)
 integrateOverIntersection computes the integral of an FE function over a grid intersection (face). More...
 
template<typename Function , class = typename Function::Space>
Function::ValueType Kaskade::integral (Function const &f)
 Evaluate integrals of finite element functions. More...
 

Detailed Description

Tools for representing and working with finite element functions.

Function Documentation

◆ forEachQP()

template<class Functor , class Space >
void Kaskade::forEachQP ( Functor &  g,
Space const &  space 
)

Loops over the whole domain occupied by the grid and calls the functor g for each integration point with an evaluator of the associated space.

g has to provide a (possibly non-const) operator()(CellPointer, QuadPoint const&, Grid::ctype, Evaluator const&) with Cell=Grid::Codim<0>::Entity, QuadPoint=QuadraturePoint<Grid::ctype,Grid::dimension> and a method int order(Evaluator const&).

WARNING: This function is deprecated in favor of an alternative interface that does not require g to have a method "order" but takes a second functor providing this order.

Definition at line 100 of file integration.hh.

Referenced by Kaskade::forEachQuadPoint(), and Kaskade::Elastomechanics::orientationPreservingStepsize().

◆ forEachQuadPoint() [1/2]

template<class Function , class Functor >
void Kaskade::forEachQuadPoint ( Function const &  f,
Functor &  accumulator 
)

Loops over the whole domain occupied by the grid and calls an accumulator for each integration point with the value of the function f. This version of forEach requires the function f to be evaluatable by an Evaluator of the associated space of f.

Parameters
fthe function to be evaluated
accumulatora functor with a (possibly non-const) operator()(Cell, QuadPoint const&, Value) with Cell=Grid::Codim<0>::Entity, QuadPoint=QuadraturePoint<Grid::ctype,Grid::dimension>, and Value=Function::ValueType.

Definition at line 73 of file integration.hh.

Referenced by Kaskade::integral(), and Kaskade::LocalIntegral< Space >::operator()().

◆ forEachQuadPoint() [2/2]

template<class Functor , class Function , class Space >
void Kaskade::forEachQuadPoint ( Function const &  f,
Functor &  g,
Space &  space 
)

Loops over the whole domain occupied by the grid and calls the functor g for each integration point with the value of the function f.

This version of forEach requires the function f to be evaluatable by cell and local coordinate.

g has to provide a (possibly non-const) operator()(Cell const&, QuadPoint const&, Value) with Cell=Grid::Codim<0>::Entity, QuadPoint=QuadraturePoint<Grid::ctype,Grid::dimension>, and Value=Function::ValueType.

Definition at line 141 of file integration.hh.

◆ integral()

template<typename Function , class = typename Function::Space>
Function::ValueType Kaskade::integral ( Function const &  f)

Evaluate integrals of finite element functions.

This computes

\[ \int_\Omega f(x)\, dx. \]

Template Parameters
afinite element function type, i.e. FunctionSpaceElement<...>, or a FunctionView.
Parameters
fa finite element function

Definition at line 289 of file integration.hh.

Referenced by Kaskade::IntegralF::id(), Kaskade::integrateOverBoundary(), Kaskade::integrateOverIntersection(), Kaskade::Integral< Space >::operator()(), Kaskade::L2Norm::square(), and Kaskade::H1SemiNorm::square().

◆ integrateOverIntersection()

template<class FEFunction >
FEFunction::ValueType Kaskade::integrateOverIntersection ( FEFunction const &  function,
typename FEFunction::Space::GridView::Intersection const &  intersection,
typename FEFunction::Space::Evaluator &  evaluator 
)

integrateOverIntersection computes the integral of an FE function over a grid intersection (face).

Parameters
functionis the FE function to be integrated.
intersectioncodim one domain of integration.
evaluatoris suitable evaluator for the function.
Template Parameters
FEFunctionis the type of the integrand.
Returns
value of integral

Definition at line 183 of file integration.hh.

Referenced by Kaskade::integrateOverBoundary().