KASKADE 7 development version
|
This evaluation cache concept defines the interface that is accessed by the assembler. More...
#include <variationalfunctional.hh>
This evaluation cache concept defines the interface that is accessed by the assembler.
This is not a base class, only a documentation!
Different objects have to be independent w.r.t. parallel execution.
Definition at line 555 of file variationalfunctional.hh.
Public Member Functions | |
DomainCache (Functional const &f, typename Vars::VariableSet const &u, int flags) | |
template<class Entity > | |
void | moveTo (Entity const &) |
template<class Position , class Evaluators > | |
void | evaluateAt (Position const &x, Evaluators const &evaluators) |
Scalar | d0 () const |
template<int row, int dim> | |
Dune::FieldVector< Scalar, Variables::template Components< row >::m > | d1 (VariationalArg< Scalar, dim > const &arg) const |
template<int row, int col, int dim> | |
Dune::FieldMatrix< Scalar, Variables::template Components< row >::m, Variables::template Components< col >::m > | d2 (VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const |
|
inline |
Constructs an evaluation cache from the associated functional f_, a point of linearization u_, and flags from the assembler (see VariationalFunctionalAssembler)
[in] | flags | A bit field that tells the cache what parts are to be assembled. If flags & Assembler::RHS is false, d1() will never be called. If flags & Assembler::VALUE is false, d0() will never be called, and d2() will never be called if flags & Assembler::MATRIX is false. This information can be used to restrict the computation to the actually required quantities. |
Definition at line 567 of file variationalfunctional.hh.
Scalar NonlinearVariationalFunctional::DomainCache::d0 | ( | ) | const |
Returns the value \(f(u(x))\)
Dune::FieldVector< Scalar, Variables::template Components< row >::m > NonlinearVariationalFunctional::DomainCache::d1 | ( | VariationalArg< Scalar, dim > const & | arg | ) | const |
Returns the directional derivative of \(f\) with respect to the variable \(u_{\mathrm{row}}\) in direction of the test function given by arg.
Dune::FieldMatrix< Scalar, Variables::template Components< row >::m, Variables::template Components< col >::m > NonlinearVariationalFunctional::DomainCache::d2 | ( | VariationalArg< Scalar, dim > const & | argT, |
VariationalArg< Scalar, dim > const & | argA | ||
) | const |
Returns the second directional derivative of \(f\) with respect to the variables \(u_{\mathrm{row}}\) and \(u_{\mathrm{col}}\) in direction of the test functions given by
This method need only be defined for values of row/col for which D2<row,col>::present is true.
|
inline |
Construct all data that is constant at one point. Default implementation in EvalCacheBase: assert(false)
Definition at line 583 of file variationalfunctional.hh.
|
inline |
Construct all data that is constant in an entity. Default implementation in EvalCacheBase: does nothing
Definition at line 574 of file variationalfunctional.hh.