KASKADE 7 development version
Public Member Functions | List of all members
VariationalFunctional::DomainCache Struct Reference

This evaluation cache concept defines the interface that is accessed by the assembler. More...

#include <variationalfunctional.hh>

Detailed Description

This evaluation cache concept defines the interface that is accessed by the assembler.

In the variational functional

\[ J(u) = \int_\Omega f(u_1(x),\nabla u_1(x),\dots,u_n(x),\nabla u_n(x)) \, dx + \int_{\partial\Omega} g(u_1(x),\dots,u_n(x)) \, dx, \]

the domain cache is responsible for evaluating the integrand \( f \) and its derivatives. The VariationalFunctionalAssembler performs the integration cell by cell. For each cell it considers, it first calls moveTo with the current cell, allowing the domain cache to conduct lengthy computations that depend only on the cell, but not on the quadrature point inside the cell (e.g., piecewise constant coefficients), and store the results for later multiple use. Next, the assembler runs through all quadrature points inside the current cell, and calls evaluateAt with the local position of the quadrature point. Again, the domain cache has the opportunity to perform expensive calculations that depend on the quadrature point, but not on the directions for which directional derivatives are to be computed (e.g., spatially non-constant coefficients of the PDE). After that, d0, d1, and d2 are called possibly several times with different variational arguments specifying the directions into which directional derivatives are to be computed.

Thus, moveTo and evaluateAt need only be implemented if there is something worthwile to compute. Otherwise, one can rely on the do-nothing default implementation of TrivialCacheBase, from which the domain cache can inherit.

This is not a base class, only a documentation!

Different objects have to be independent w.r.t. parallel execution.

Definition at line 142 of file variationalfunctional.hh.

Public Member Functions

void moveTo (Cell const &)
 Construct all data that is constant in an entity. Default implementation in EvalCacheBase: does nothing. More...
 
template<class Position , class Evaluators >
void evaluateAt (Position const &x, Evaluators const &evaluators)
 Construct all data that is constant at one point. Default implementation in EvalCacheBase: assert(false) More...
 
Scalar d0 () const
 
template<int row, int dim>
Dune::FieldVector< Scalar, TestVars::template Components< row >::m > d1 (VariationalArg< Scalar, dim > const &arg) const
 
template<int row, int col, int dim>
Dune::FieldMatrix< Scalar, TestVars::template Components< row >::m, AnsatzVars::template Components< row >::m > d2 (VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
 

Member Function Documentation

◆ d0()

Scalar VariationalFunctional::DomainCache::d0 ( ) const

Returns the value \(f(0)\)

◆ d1()

template<int row, int dim>
Dune::FieldVector< Scalar, TestVars::template Components< row >::m > VariationalFunctional::DomainCache::d1 ( VariationalArg< Scalar, dim > const &  arg) const

Returns the directional derivative of \(f\) at \( 0 \) with respect to the variable \( u_{\mathrm{row}}\) in direction of the test function given by arg.

◆ d2()

template<int row, int col, int dim>
Dune::FieldMatrix< Scalar, TestVars::template Components< row >::m, AnsatzVars::template Components< row >::m > VariationalFunctional::DomainCache::d2 ( VariationalArg< Scalar, dim > const &  argT,
VariationalArg< Scalar, dim > const &  argA 
) const

Returns the second directional derivative of \(f\) at \( 0 \) with respect to the variables \(u_{\mathrm{row}}\) and \(u_{\mathrm{col}}\) in direction of the test functions given by argT and argA.

For multi-component variables (e.g., displacement in solid mechanics), the given test and ansatz variations argT and argA will usually be scalar. The returned matrix \(R\) is then defined entry-wise as

\[ R_{ij} = f''(0)[\mathrm{argT}e_i,\mathrm{argA}e_j], \]

i.e. every return entry corresponds to a second directional derivative in directions where just one component is changed in each.

This method need only be defined for values of row/col for which D2<row,col>::present is true.

◆ evaluateAt()

template<class Position , class Evaluators >
void VariationalFunctional::DomainCache::evaluateAt ( Position const &  x,
Evaluators const &  evaluators 
)
inline

Construct all data that is constant at one point. Default implementation in EvalCacheBase: assert(false)

Definition at line 158 of file variationalfunctional.hh.

◆ moveTo()

void VariationalFunctional::DomainCache::moveTo ( Cell const &  )
inline

Construct all data that is constant in an entity. Default implementation in EvalCacheBase: does nothing.

Definition at line 149 of file variationalfunctional.hh.


The documentation for this struct was generated from the following file: