KASKADE 7 development version
Public Types | Public Member Functions | List of all members
NonlinearVariationalFunctional::BoundaryCache Struct Reference

Defining boundary conditions. More...

#include <variationalfunctional.hh>

Detailed Description

Defining boundary conditions.

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 624 of file variationalfunctional.hh.

Inheritance diagram for NonlinearVariationalFunctional::BoundaryCache:

Public Types

typedef AnsatzVars::Grid::LeafIntersectionIterator FaceIterator
 Type of the Face Information. More...
 

Public Member Functions

 BoundaryCache (Functional const &f, typename Vars::VariableSet const &u, int const flags=7)
 
template<class FaceIterator >
void moveTo (FaceIterator 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
 

Member Typedef Documentation

◆ FaceIterator

typedef AnsatzVars::Grid::LeafIntersectionIterator NonlinearVariationalFunctional::BoundaryCache::FaceIterator

Type of the Face Information.

Definition at line 629 of file variationalfunctional.hh.

Constructor & Destructor Documentation

◆ BoundaryCache()

NonlinearVariationalFunctional::BoundaryCache::BoundaryCache ( Functional const &  f,
typename Vars::VariableSet const &  u,
int const  flags = 7 
)
inline

Constructs an evaluation cache from the associated functional f_, a point of linearization u_, and flags from the assembler

Parameters
[in]flagsA 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 640 of file variationalfunctional.hh.

Member Function Documentation

◆ d0()

Scalar NonlinearVariationalFunctional::BoundaryCache::d0 ( ) const

Returns the value \(g(u(x))\)

◆ d1()

template<int row, int dim>
Dune::FieldVector< Scalar, Variables::template Components< row >::m > NonlinearVariationalFunctional::BoundaryCache::d1 ( VariationalArg< Scalar, dim > const &  arg) const

Returns the directional derivative of \(g\), 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, Variables::template Components< row >::m, Variables::template Components< col >::m > NonlinearVariationalFunctional::BoundaryCache::d2 ( VariationalArg< Scalar, dim > const &  argT,
VariationalArg< Scalar, dim > const &  argA 
) const

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

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 NonlinearVariationalFunctional::BoundaryCache::evaluateAt ( Position const &  x,
Evaluators const &  evaluators 
)
inline

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

Definition at line 654 of file variationalfunctional.hh.

◆ moveTo()

template<class FaceIterator >
void NonlinearVariationalFunctional::BoundaryCache::moveTo ( FaceIterator const &  )
inline

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

Definition at line 646 of file variationalfunctional.hh.


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