KASKADE 7 development version
Classes | Public Types | Public Member Functions | Static Public Attributes | Protected Attributes | Related Functions | List of all members
Kaskade::LinearizationAt< Functional_ > Class Template Reference

Proxy class for the linearization of a nonlinear functional. More...

#include <functional_aux.hh>

Detailed Description

template<class Functional_>
class Kaskade::LinearizationAt< Functional_ >

Proxy class for the linearization of a nonlinear functional.

This adaptor maps the NonlinearVariationalFunctional concept to the VariationalFunctional concept.

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 a variable \( u\in U_o \) one has to create a linearization, which is represented by this class.

Definition at line 614 of file functional_aux.hh.

Inheritance diagram for Kaskade::LinearizationAt< Functional_ >:

Classes

struct  D1
 Rhs block information. More...
 
struct  D2
 matrix block information More...
 

Public Types

typedef Functional_ Functional
 
typedef Functional::Scalar Scalar
 
typedef Functional::AnsatzVars AnsatzVars
 
typedef Functional::TestVars TestVars
 
typedef Functional::OriginVars OriginVars
 
typedef Functional::DomainCache DomainCache
 
typedef Functional::BoundaryCache BoundaryCache
 
typedef Functional::OriginVars::VariableSet PointOfLinearization
 Type of variables around which the functional will be linearized. More...
 

Public Member Functions

 LinearizationAt (Functional const &f_, PointOfLinearization const &u_)
 Construct a linearization 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...
 
template<int row, int col>
auto makePositiveThreshold () const
 
PointOfLinearization const & getOrigin () const
 Get point of linearization. More...
 
Functional const & getFunctional () const
 Get functional. More...
 
int derivatives () const
 Gives the highest derivative that arises in the weak formulation. More...
 

Static Public Attributes

static ProblemType const type = Functional::type
 

Protected Attributes

Functional const & f
 
PointOfLinearization const & u
 

Related Functions

(Note that these are not member functions.)

template<class Functional >
LinearizationAt< Functionallinearization (Functional const &f, typename Functional::OriginVars::VariableSet const &u)
 Convenience routine: construct linearization without having to know the type of the functional. More...
 

Member Typedef Documentation

◆ AnsatzVars

template<class Functional_ >
typedef Functional::AnsatzVars Kaskade::LinearizationAt< Functional_ >::AnsatzVars

Definition at line 622 of file functional_aux.hh.

◆ BoundaryCache

template<class Functional_ >
typedef Functional::BoundaryCache Kaskade::LinearizationAt< Functional_ >::BoundaryCache

Definition at line 627 of file functional_aux.hh.

◆ DomainCache

template<class Functional_ >
typedef Functional::DomainCache Kaskade::LinearizationAt< Functional_ >::DomainCache

Definition at line 626 of file functional_aux.hh.

◆ Functional

template<class Functional_ >
typedef Functional_ Kaskade::LinearizationAt< Functional_ >::Functional

Definition at line 617 of file functional_aux.hh.

◆ OriginVars

template<class Functional_ >
typedef Functional::OriginVars Kaskade::LinearizationAt< Functional_ >::OriginVars

Definition at line 624 of file functional_aux.hh.

◆ PointOfLinearization

template<class Functional_ >
typedef Functional::OriginVars::VariableSet Kaskade::LinearizationAt< 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. In any case, it must refer to the same list of FE spaces.

\TODO: Why that? 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 641 of file functional_aux.hh.

◆ Scalar

template<class Functional_ >
typedef Functional::Scalar Kaskade::LinearizationAt< Functional_ >::Scalar

Definition at line 620 of file functional_aux.hh.

◆ TestVars

template<class Functional_ >
typedef Functional::TestVars Kaskade::LinearizationAt< Functional_ >::TestVars

Definition at line 623 of file functional_aux.hh.

Constructor & Destructor Documentation

◆ LinearizationAt()

template<class Functional_ >
Kaskade::LinearizationAt< Functional_ >::LinearizationAt ( Functional const &  f_,
PointOfLinearization const &  u_ 
)
inline

Construct a linearization from a given functional.

Definition at line 646 of file functional_aux.hh.

Member Function Documentation

◆ createBoundaryCache()

template<class Functional_ >
BoundaryCache Kaskade::LinearizationAt< Functional_ >::createBoundaryCache ( int  flags) const
inline

create a boundary cache

Definition at line 658 of file functional_aux.hh.

◆ createDomainCache()

template<class Functional_ >
DomainCache Kaskade::LinearizationAt< Functional_ >::createDomainCache ( int  flags) const
inline

create a domain cache

Definition at line 652 of file functional_aux.hh.

◆ derivatives()

template<class Functional_ >
int Kaskade::LinearizationAt< Functional_ >::derivatives ( ) const
inline

Gives the highest derivative that arises in the weak formulation.

Definition at line 703 of file functional_aux.hh.

◆ getFunctional()

template<class Functional_ >
Functional const & Kaskade::LinearizationAt< Functional_ >::getFunctional ( ) const
inline

Get functional.

Definition at line 695 of file functional_aux.hh.

◆ getOrigin()

template<class Functional_ >
PointOfLinearization const & Kaskade::LinearizationAt< Functional_ >::getOrigin ( ) const
inline

Get point of linearization.

Definition at line 689 of file functional_aux.hh.

◆ integrationOrder()

template<class Functional_ >
template<class Cell >
int Kaskade::LinearizationAt< Functional_ >::integrationOrder ( Cell const &  cell,
int  shapeFunctionOrder,
bool  boundary 
) const
inline

integration order

Definition at line 677 of file functional_aux.hh.

◆ makePositiveThreshold()

template<class Functional_ >
template<int row, int col>
auto Kaskade::LinearizationAt< Functional_ >::makePositiveThreshold ( ) const
inline

Definition at line 683 of file functional_aux.hh.

Member Data Documentation

◆ f

template<class Functional_ >
Functional const& Kaskade::LinearizationAt< Functional_ >::f
protected

◆ type

template<class Functional_ >
ProblemType const Kaskade::LinearizationAt< Functional_ >::type = Functional::type
static

Definition at line 618 of file functional_aux.hh.

◆ u

template<class Functional_ >
PointOfLinearization const& Kaskade::LinearizationAt< Functional_ >::u
protected

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