KASKADE 7 development version
Classes | Public Types | Public Member Functions | Static Public Attributes | List of all members
VariationalFunctional Class Reference

Documentation of the concept of a quadratic variational functionalThe variational functional concept defines the interface that is accessed by the VariationalFunctionalAssembler. More...

#include <variationalfunctional.hh>

Detailed Description

Documentation of the concept of a quadratic variational functional

The variational functional concept defines the interface that is accessed by the VariationalFunctionalAssembler.

This is not a base class, only a documentation! A convenience base class defining more or less useful default values for several required member types and functions is FunctionalBase.

The mathematical concept represented by this interface is that of a 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)) \, ds + \sum_{F_T\in\mathcal{F}} g_F(u_1(x^+),\dots,u_n(x^+),u_1(x^-),\dots,u_n(x^-)) \, ds \]

with quadratic functions \(f\), \(g\), and \(g_F\) to be evaluated at \(0\). The different variables \(u_i\) may be scalar or vector-valued, and may live in different finite element function spaces. \( \mathcal{F} \) is a (possibly empty) subset of cell-face incidences in the grid, and \(x^+\) and \(x^-\) refer to the same point on different sides of the face. This way, jumps of discontinuous ansatz functions can be integrated on faces.

Some local classes have to be written into VariationalFunctional, where most of the functionality of the functional is defined. These are:

Also read the documentation of LinearizationAt in functional_aux.hh, an adapter class for the usage and the definition of a variational functional, and the NonlinearVariationalFunctional concept. If only quadratic functionals (linear problems) are to be defined, these classes are not necessary.

Const methods are required to be thread-safe. Note that the space lists AnsatzVars::Spaces and TestVars::Spaces must coincide!

See also
LinearizationAt
NonlinearVariationalFunctional

Definition at line 63 of file variationalfunctional.hh.

Classes

struct  BoundaryCache
 Evaluation of boundary conditions. More...
 
struct  D1
 Provides static information about the right hand side blocks. More...
 
struct  D2
 Provides static information about the submatrix blocks of the stiffness block matrix. More...
 
struct  DomainCache
 This evaluation cache concept defines the interface that is accessed by the assembler. More...
 
struct  InnerBoundaryCache
 Evaluates jump contributions at interior faces, suitable for discontinuous Galerkin methods (optional). More...
 

Public Types

typedef unspecified Scalar
 The scalar type to be used for this variational problem. More...
 
typedef unspecified AnsatzVars
 A description of the ansatz variables occuring in the variational problem. More...
 
typedef unspecified TestVars
 A description of the test variables occuring in the "variational" problem. More...
 
typedef unspecified OriginVars
 A description of the variables defining the linearization point. More...
 
typedef Vars::Grid::template Codim< 0 >::Entity Entity
 

Public Member Functions

DomainCache createDomainCache (int flags) const
 Creates a DomainCache. More...
 
BoundaryCache createBoundaryCache (int flags) const
 Creates a BoundaryCache. More...
 
InnerBoundaryCache createInnerBoundaryCache (int flags) const
 Creates an InnerBoundaryCache. More...
 
int integrationOrder (typename Variables::Grid::template Codim< 0 >::Entity const &cell, int shapeFunctionOrder, bool boundary) const
 This method defines a suitable quadrature order to be used on the given cell. More...
 
bool considerFace (Face< typename Variables::GridView > const &face) const
 

Static Public Attributes

static ProblemType const type
 The type of problem, either VariationalFunctional or WeakFormulation. More...
 

Member Typedef Documentation

◆ AnsatzVars

typedef unspecified VariationalFunctional::AnsatzVars

A description of the ansatz variables occuring in the variational problem.

This should be an instance of VariableSetDescription<...>.

Definition at line 76 of file variationalfunctional.hh.

◆ Entity

typedef Vars::Grid::template Codim<0>::Entity VariationalFunctional::Entity

The following typedef should be useful, when accessing Entities, it is however not strictly required for a functional to work

Definition at line 110 of file variationalfunctional.hh.

◆ OriginVars

typedef unspecified VariationalFunctional::OriginVars

A description of the variables defining the linearization point.

In most cases, this is simply AnsatzVars.

In some cases, however, the linear(ized) problem to be solved is posed in a different ansatz space than the current evaluation point, e.g. hierarchical error estimators computing corrections in an extension space. Then there is a difference between ansatz and test variables on one hand and the origin variables on the other.

This should be an instance of VariableSetDescription<...>.

Definition at line 99 of file variationalfunctional.hh.

◆ Scalar

typedef unspecified VariationalFunctional::Scalar

The scalar type to be used for this variational problem.

Definition at line 69 of file variationalfunctional.hh.

◆ TestVars

typedef unspecified VariationalFunctional::TestVars

A description of the test variables occuring in the "variational" problem.

This should be an instance of VariableSetDescription<...>.

Definition at line 84 of file variationalfunctional.hh.

Member Function Documentation

◆ considerFace()

bool VariationalFunctional::considerFace ( Face< typename Variables::GridView > const &  face) const

This method tells on which inner cell-face incidences \( F_T \) of the grid view the integration of the variational functional shall be performed, if inner faces shall be treated at all, i.e. an InnerBoundaryCache class is present. The stiffness matrix entries are only allocated for those faces that are reported true, and the assembly is restricted to those cells. Thus, a safe but maybe inefficient version is to always return true, which is the default implementation in FunctionalBase.

The method needs to be implemented only if the class InnerBoundaryCache is defined.

◆ createBoundaryCache()

BoundaryCache VariationalFunctional::createBoundaryCache ( int  flags) const
inline

Creates a BoundaryCache.

Parameters
[in]flagsa bit field describing what values will be accessed. See Assembler enums. This should be passed on to the boundary cache, which then can restrict the computation to values that are actually needed.

Definition at line 429 of file variationalfunctional.hh.

◆ createDomainCache()

DomainCache VariationalFunctional::createDomainCache ( int  flags) const
inline

Creates a DomainCache.

Parameters
[in]flagsa bit field describing what values will be accessed. See Assembler enums. This should be passed on to the domain cache, which then can restrict the computation to values that are actually needed.

Definition at line 420 of file variationalfunctional.hh.

◆ createInnerBoundaryCache()

InnerBoundaryCache VariationalFunctional::createInnerBoundaryCache ( int  flags) const
inline

Creates an InnerBoundaryCache.

Parameters
[in]flagsa bit field describing what values will be accessed. See Assembler enums. This should be passed on to the inner boundary cache, which then can restrict the computation to values that are actually needed.

This is optional and has to be implemented only if the InnerBoundaryCache structure is defined.

Definition at line 440 of file variationalfunctional.hh.

◆ integrationOrder()

int VariationalFunctional::integrationOrder ( typename Variables::Grid::template Codim< 0 >::Entity const &  cell,
int  shapeFunctionOrder,
bool  boundary 
) const

This method defines a suitable quadrature order to be used on the given cell.

The maximum shape function order occuring in the finite element spaces for this cell is provided for convenience. If boundary is true, this refers to evaluation of \(g\), otherwise to the evaluation of \(f\).

Member Data Documentation

◆ type

ProblemType const VariationalFunctional::type
static

The type of problem, either VariationalFunctional or WeakFormulation.

Definition at line 104 of file variationalfunctional.hh.


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