KASKADE 7 development version
|
Proxy class for evaluating differences of linearizations of a nonlinear functional. More...
#include <functional_aux.hh>
Proxy class for evaluating differences of linearizations of a nonlinear functional.
Acceptance tests in Newton's method for minimization problems usually involves some kind of comparison \( f(u+\delta u) < f(u) \), where
\[ f(u) = \sum_{i=1}^N f_i(u). \]
Now taking the difference of both sums (large, because it's assembly over the grid) is numerically unstable for \( \delta u \) small. A more stable evaluation is
\[ \sum_{i=1}^N (f_i(u)-f_i(u+\delta u)) > 0. \]
This class supports the evaluation of the sum of differences.
According to the following idea: A VariationalFunctional defines a mapping \( F : U \to V \), which can be evaluated and linearized in an appropriate ansatz space \( U_a\subset U \), only if an argument \( u \in U_o \) in a - possibly different - "origin" subspace \( U_o \subset U \) is given.
Hence, from a functional \( f \) and two variables \( u_1, u_2 \in U_o \) one has to create a linearization of \( f(u_1)-f(u_2) \), which is represented by this class.
Definition at line 748 of file functional_aux.hh.
Classes | |
class | BoundaryCache |
struct | D1 |
Rhs block information. More... | |
struct | D2 |
matrix block information More... | |
class | DomainCache |
Public Types | |
typedef Functional_ | Functional |
typedef Functional::Scalar | Scalar |
typedef Functional::AnsatzVars | AnsatzVars |
typedef Functional::TestVars | TestVars |
typedef Functional::OriginVars | OriginVars |
typedef Functional::OriginVars::VariableSet | PointOfLinearization |
Type of variables around which the functional will be linearized. More... | |
Public Member Functions | |
LinearizationDifferenceAt (Functional const &f_, PointOfLinearization const &u1_, PointOfLinearization const &u2_) | |
Construct a difference of linearizations from a given functional. More... | |
DomainCache | createDomainCache (int flags) const |
create a domain cache More... | |
BoundaryCache | createBoundaryCache (int flags) const |
create a boundary cache More... | |
template<class Cell > | |
int | integrationOrder (Cell const &cell, int shapeFunctionOrder, bool boundary) const |
integration order More... | |
Functional const & | getFunctional () const |
Get functional. More... | |
Static Public Attributes | |
static ProblemType const | type = Functional::type |
Protected Attributes | |
Functional const & | f |
PointOfLinearization const & | u1 |
PointOfLinearization const & | u2 |
typedef Functional::AnsatzVars Kaskade::LinearizationDifferenceAt< Functional_ >::AnsatzVars |
Definition at line 756 of file functional_aux.hh.
typedef Functional_ Kaskade::LinearizationDifferenceAt< Functional_ >::Functional |
Definition at line 751 of file functional_aux.hh.
typedef Functional::OriginVars Kaskade::LinearizationDifferenceAt< Functional_ >::OriginVars |
Definition at line 758 of file functional_aux.hh.
typedef Functional::OriginVars::VariableSet Kaskade::LinearizationDifferenceAt< Functional_ >::PointOfLinearization |
Type of variables around which the functional will be linearized.
The type of "iterate", i.e. the point around which the functional is linearized in the ansatz and test spaces. Note that this can be, but need not be, the ansatz space itself. It must have the same number of variables with the same number of components each, but the underlying FE spaces can be different.
Definition at line 770 of file functional_aux.hh.
typedef Functional::Scalar Kaskade::LinearizationDifferenceAt< Functional_ >::Scalar |
Definition at line 754 of file functional_aux.hh.
typedef Functional::TestVars Kaskade::LinearizationDifferenceAt< Functional_ >::TestVars |
Definition at line 757 of file functional_aux.hh.
|
inline |
Construct a difference of linearizations from a given functional.
Definition at line 775 of file functional_aux.hh.
|
inline |
create a boundary cache
Definition at line 850 of file functional_aux.hh.
|
inline |
create a domain cache
Definition at line 844 of file functional_aux.hh.
|
inline |
Get functional.
Definition at line 873 of file functional_aux.hh.
|
inline |
integration order
Definition at line 867 of file functional_aux.hh.
|
protected |
Definition at line 876 of file functional_aux.hh.
Referenced by Kaskade::LinearizationDifferenceAt< Functional_ >::createBoundaryCache(), Kaskade::LinearizationDifferenceAt< Functional_ >::createDomainCache(), Kaskade::LinearizationDifferenceAt< Functional_ >::getFunctional(), and Kaskade::LinearizationDifferenceAt< Functional_ >::integrationOrder().
|
static |
Definition at line 752 of file functional_aux.hh.
|
protected |
Definition at line 877 of file functional_aux.hh.
Referenced by Kaskade::LinearizationDifferenceAt< Functional_ >::createBoundaryCache(), and Kaskade::LinearizationDifferenceAt< Functional_ >::createDomainCache().
|
protected |
Definition at line 878 of file functional_aux.hh.
Referenced by Kaskade::LinearizationDifferenceAt< Functional_ >::createBoundaryCache(), and Kaskade::LinearizationDifferenceAt< Functional_ >::createDomainCache().