KASKADE 7 development version
|
Evaluates jump contributions at interior faces, suitable for discontinuous Galerkin methods (optional). More...
#include <variationalfunctional.hh>
Evaluates jump contributions at interior faces, suitable for discontinuous Galerkin methods (optional).
If this structure is defined, the assembler does assemble jump terms
\[ \int_F g(u^{\rm left},u^{\rm right})\,ds \]
at interior faces \( F \) of the grid. If no InnerBoundaryCache class is defined, the assembler ignores these integrals.
Definition at line 259 of file variationalfunctional.hh.
Public Member Functions | |
template<class FaceIterator > | |
void | moveTo (FaceIterator const &) |
template<class Position , class Evaluators > | |
void | evaluateAt (Position const &x, Evaluators const &evaluators, Evaluators const &neighbourEvaluators) |
Prepare the cache for evaluation at given quadrature point. More... | |
Scalar | d0 () const |
Returns the value \(g(0)\). More... | |
template<int row, int dim> | |
Dune::FieldVector< Scalar, TestVars::template Components< row >::m > | d1 (VariationalArg< Scalar, dim > const &argT) 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, bool centerCell) const |
Scalar VariationalFunctional::InnerBoundaryCache::d0 | ( | ) | const |
Returns the value \(g(0)\).
Dune::FieldVector< Scalar, TestVars::template Components< row >::m > VariationalFunctional::InnerBoundaryCache::d1 | ( | VariationalArg< Scalar, dim > const & | argT | ) | const |
Returns the directional derivative of \(g\) at \( 0 \), with respect to the variable \(u_{\mathrm{row}}\) in direction of the test function given by arg.
The test variational argument belongs to the current cell, i.e. the interior side of the face.
Dune::FieldMatrix< Scalar, TestVars::::template Components< row >::m, AnsatzVars::::template Components< row >::m > VariationalFunctional::InnerBoundaryCache::d2 | ( | VariationalArg< Scalar, dim > const & | argT, |
VariationalArg< Scalar, dim > const & | argA, | ||
bool | centerCell | ||
) | const |
Returns the second directional derivative of \( g \) at \( 0 \), 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.
If centerCell is true, both test and ansatz variational arguments are associated on the current cell. Otherwise, the ansatz argument argA belongs to the neighbour cell, i.e. the outside of the face.
This method need only be defined for values of row/col for which D2<row,col>::present is true.
void VariationalFunctional::InnerBoundaryCache::evaluateAt | ( | Position const & | x, |
Evaluators const & | evaluators, | ||
Evaluators const & | neighbourEvaluators | ||
) |
Prepare the cache for evaluation at given quadrature point.
x | the quadrature point in the local face coordinate system |
evaluators | evaluators for quantities on the current cell |
neighbourEvaluators | evaluators for the cell on the opposite side of the face |
void VariationalFunctional::InnerBoundaryCache::moveTo | ( | FaceIterator const & | ) |