KASKADE 7 development version
|
Linearly implicit Euler method. More...
#include <semieuler.hh>
Linearly implicit Euler method.
A class that defines the weak formulation of the stationary elliptic problem that results from the linearly implicit Euler method for
\[ B \dot u = L(u) , \]
which is, doing the linearization at \( u_0 \),
\[ (B - \tau L'(u_0)) \delta u_k = \tau L(u_k), \quad u_{k+1} = u_k + \delta u_k. \]
The left hand side \( B - \tau L'(u_0) \) is assembled as matrix, and the right hand side \( \tau L(u_k) \) is assembled as right hand side by the assembler, if fed with this weak problem formulation:
\( B \) may depend on \( u \), made explicit by setting
in the parabolic equation class. In this case, \( B(u) \) has to be an invertable Nemyckii operator (a diagonal operator, mapping \( u(x) \) to \( B(x,u(x)) u(x) \). Then the linearly implicit Euler scheme is (see below for derivation)
\[ (B(u_0) - \tau L'(u_0) + \tau_k B^{-1}(u_0)B'(u_0)\mathrm{diag}L(u_0)) \delta u_k = \tau_k L(u)+\frac{\tau_k}{\tau_{k-1}}(B(u_0)-B(u_k))\delta u_{k-1}. \]
The values \( u_0 \) (linearization point for the Jacobian), \( u_k \) (evaluation point for the right hand side), and \( \tau_k/\tau_{k-1}\, \delta u_{k-1} \) (previous increment) have to be specified when constructing the DomainCache. This is most conveniently done via SemiLinearizationAt.
PE | the weak formulation of the parabolic equation, adhering to the ParabolicEquation interface. |
Derivation: Lubich and Roche (Rosenbrock Methods for Differential-algebraic Systems with Solution-dependent Singular Matrix Multiplying the Derivative, Computing 43:325-342, 1990) suggest to reformulate the parabolic system with non-constant \( B \) above as
\[ \begin{aligned} \dot u &= z \\ 0 &= L(u)-B(u)z \end{aligned} \]
and applying, e.g., Rosenbrock methods. Here, simply a linearly implicit Euler. Linearizing the right hand side of the extended system at \( u_0, z_0 \) yields
\[ \Bigg( \begin{bmatrix} I & 0 \\ 0 & 0 \end{bmatrix} - \tau_k \begin{bmatrix} 0 & I \\ L'(u_0)-B'(u_0)z_0 & -B(u_0) \end{bmatrix} \Bigg) \begin{bmatrix} \delta u_k \\ \delta z_k \end{bmatrix} = \tau_k \begin{bmatrix} z_k \\ L(u_k)-B(u_k)z_k \end{bmatrix}. \]
The first row immediately gives \( \delta u_k = \tau_k (z_k + \delta z_k) \), and the second one
\[ - \tau_k (L'(u_0)-B'(u_0)z_0) \delta u_k + \tau (B(u_0)\delta z_k + B(u_k)z_k) = \tau_k f(u_k). \]
The second term on the left hand side is just
\[ \tau_k (B(u_0)(z_k+\delta z_k) + (B(u_k)-B(u_0))z_k) = B(u_0)\delta u_k + (B(u_k)-B(u_0)) \tau_k z_k, \]
hence we obtain
\[ \left(B(u_0) - \tau_k(L'(u_0)-B'(u_0)z_0)\right)\delta u_k = \tau_k L(u_k)+(B(u_0)-B(u_k)) \frac{\tau_k}{\tau_{k-1}}\delta u_{k-1}. \]
Choosing a consistent initial value \( z_0 = B(u_0)^{-1} L(u_0) \) finally results in the form given above.
Definition at line 197 of file semieuler.hh.
Classes | |
class | BoundaryCache |
struct | D1 |
struct | D2 |
class | DomainCache |
Public Types | |
typedef PE | ParabolicEquation |
typedef ParabolicEquation::Scalar | Scalar |
typedef ParabolicEquation::AnsatzVars | AnsatzVars |
typedef ParabolicEquation::TestVars | TestVars |
typedef ParabolicEquation::OriginVars | OriginVars |
typedef AnsatzVars::Grid | Grid |
using | GridView = typename AnsatzVars::GridView |
typedef PE | EvolutionEquation |
typedef EvolutionEquation::AnsatzVars::template CoefficientVectorRepresentation< 0, EvolutionEquation::TestVars::noOfVariables >::type | CoefficientVectors |
Public Member Functions | |
SemiImplicitEulerStep (ParabolicEquation *eq_, Scalar dt) | |
Constructor. More... | |
template<class Cell > | |
int | integrationOrder (Cell const &cell, int shapeFunctionOrder, bool boundary) const |
Self & | setTau (Scalar dt) |
Sets the time step size. More... | |
Scalar | getTau () const |
Access to the underlying parabolic equation. | |
ParabolicEquation & | parabolicEquation () |
ParabolicEquation const & | parabolicEquation () const |
Static Public Attributes | |
static ProblemType const | type = WeakFormulation |
typedef ParabolicEquation::AnsatzVars Kaskade::SemiImplicitEulerStep< PE >::AnsatzVars |
Definition at line 208 of file semieuler.hh.
typedef EvolutionEquation::AnsatzVars::template CoefficientVectorRepresentation<0,EvolutionEquation::TestVars::noOfVariables>::type Kaskade::SemiImplicitEulerStep< PE >::CoefficientVectors |
Definition at line 218 of file semieuler.hh.
typedef PE Kaskade::SemiImplicitEulerStep< PE >::EvolutionEquation |
Definition at line 216 of file semieuler.hh.
typedef AnsatzVars::Grid Kaskade::SemiImplicitEulerStep< PE >::Grid |
Definition at line 211 of file semieuler.hh.
using Kaskade::SemiImplicitEulerStep< PE >::GridView = typename AnsatzVars::GridView |
Definition at line 212 of file semieuler.hh.
typedef ParabolicEquation::OriginVars Kaskade::SemiImplicitEulerStep< PE >::OriginVars |
Definition at line 210 of file semieuler.hh.
typedef PE Kaskade::SemiImplicitEulerStep< PE >::ParabolicEquation |
Definition at line 205 of file semieuler.hh.
typedef ParabolicEquation::Scalar Kaskade::SemiImplicitEulerStep< PE >::Scalar |
Definition at line 207 of file semieuler.hh.
typedef ParabolicEquation::TestVars Kaskade::SemiImplicitEulerStep< PE >::TestVars |
Definition at line 209 of file semieuler.hh.
|
inline |
Constructor.
eq | a pointer to a parabolic equation object. This needs to exist during the life time of the SemiImplicitEulerStep object. |
dt | the time step. dt>=0 has to hold. |
Definition at line 412 of file semieuler.hh.
|
inline |
|
inline |
Definition at line 418 of file semieuler.hh.
|
inline |
Definition at line 442 of file semieuler.hh.
Referenced by Kaskade::Limex< Eq >::advanceTime(), Kaskade::Limex< Eq >::estimateError(), and Kaskade::Limex< Eq >::step().
|
inline |
Definition at line 443 of file semieuler.hh.
|
inline |
Sets the time step size.
dt | The new time step size. dt>=0 has to hold. |
Definition at line 428 of file semieuler.hh.
Referenced by Kaskade::Limex< Eq >::step().
|
static |
Definition at line 214 of file semieuler.hh.