KASKADE 7 development version
|
Documentation of the concept of a nonlinear variational functional The nonlinear variational functional concept defines the interface that is accessed by the LinearizationAt and linearizationAt adapters. It differs from the basic VariationalFunctional concept mainly in how the DomainCache and the BoundaryCache are constructed. More...
#include <variationalfunctional.hh>
Documentation of the concept of a nonlinear variational functional The nonlinear variational functional concept defines the interface that is accessed by the LinearizationAt and linearizationAt adapters. It differs from the basic VariationalFunctional concept mainly in how the DomainCache and the BoundaryCache are constructed.
This is not a base class, only a documentation!
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)) \, dx \]
to be evaluated for a specific fixed value of \(u\). The different variables \(u_i\) may be scalar or vector-valued, and may live in different finite element function spaces.
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, which provides helper classes for the usage and the definition of a variational functional. 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!
Definition at line 508 of file variationalfunctional.hh.
Classes | |
struct | BoundaryCache |
Defining boundary conditions. More... | |
struct | D1 |
struct | D2 |
struct | DomainCache |
This evaluation cache concept defines the interface that is accessed by the assembler. 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. This should be an instance of VariableSetDescription<...>. More... | |
typedef unspecified | TestVars |
typedef unspecified | OriginVars |
A description of the type of the point of linearization. More... | |
typedef Vars::Grid::template Codim< 0 >::Entity | Entity |
Public Member Functions | |
int | integrationOrder (typename Variables::Grid::template Codim< 0 >::Entity const &cell, int shapeFunctionOrder, bool boundary) const |
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... | |
typedef unspecified NonlinearVariationalFunctional::AnsatzVars |
A description of the ansatz variables occuring in the variational problem. This should be an instance of VariableSetDescription<...>.
Definition at line 521 of file variationalfunctional.hh.
typedef Vars::Grid::template Codim<0>::Entity NonlinearVariationalFunctional::Entity |
The following typedef should be useful, when accessing Entities, it is however not strictly required for a functional to work
Definition at line 546 of file variationalfunctional.hh.
typedef unspecified NonlinearVariationalFunctional::OriginVars |
A description of the type of the point of linearization.
Definition at line 534 of file variationalfunctional.hh.
typedef unspecified NonlinearVariationalFunctional::Scalar |
The scalar type to be used for this variational problem.
Definition at line 514 of file variationalfunctional.hh.
typedef unspecified NonlinearVariationalFunctional::TestVars |
A description of the test variables occuring in the "variational" problem. This should be an instance of VariableSetDescription<...>.
Definition at line 529 of file variationalfunctional.hh.
bool NonlinearVariationalFunctional::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.
int NonlinearVariationalFunctional::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\).
|
static |
The type of problem, either VariationalFunctional or WeakFormulation.
Definition at line 540 of file variationalfunctional.hh.