KASKADE 7 development version
Public Member Functions | List of all members
Kaskade::Sdc< Vector > Class Template Reference

Encapsulates the state of an SDC fixed point iteration and provides stepping. More...

#include <sdc.hh>

Detailed Description

template<class Vector>
class Kaskade::Sdc< Vector >

Encapsulates the state of an SDC fixed point iteration and provides stepping.

This solves the initial value problem \( \dot u = f(u) \) on the interval \( [t_0,T] \) by spectral deferred correction. The equivalent Picard equation

\[ u(t) = u(a) + \int_{\tau=t_0}^t f(u(\tau)) \, d\tau \]

is considered, and approximated by discretizing \( u\in S\times[t_0,T] \) by some polynomial \( P_{\le N}[S] \) and replacing the integral by quadrature on collocation points \( t_1, \dots, t_N \in \mathopen]t_0,T]\). Representing the polynomial \( u \) by its values \( u_i \) at \( t_i \), \( i=0,\dots, N\), we obtain

\[ u_{i+1} = u_i + \sum_j S_{ij} f(u_j). \]

Given some approximation \( \hat u \), the error \( \delta u = u - \hat u\) satisfies

\[ \begin{aligned} \delta u_{i+1} &= u_i + \sum_j S_{ij} f(u_j) - \hat u_i + (\hat u_i - \hat u_{i+1}) \\ &= \delta u_i + \sum_j S_{ij} f(\hat u_j + \delta u_j) + (\hat u_i - \hat u_{i+1}) \end{aligned}. \]

Linearizing \( f \) around \( \hat u_j \) and splitting the sum leads to

\[ \delta u_{i+1} = \delta u_i + \sum_j S_{ij} f(\hat u_j) - (\hat u_{j+1} - \hat u_j) + \sum_j S_{ij} f'(\hat u_j) \delta u_j, \]

which can be simplified by replacing \( S \) by a lower triangular approximation \( \hat S \). This defines a correction \( \hat \delta u \) and consequently a fixed point iteration:

\[ \begin{gathered}\hat \delta u_{i+1} := \hat\delta u_i + \sum_j S_{ij} f(\hat u_j) - (\hat u_{j+1} - \hat u_j) + \sum_{j\le i+1} \hat S_{ij} f'(\hat u_j) \hat\delta u_j, \\ \hat u \gets \hat u + \hat \delta u. \end{gathered} \]

Template Parameters
Vectortype of vectors from some linear space.

Definition at line 1059 of file sdc.hh.

Public Member Functions

 Sdc (SDCTimeGrid const &ts_, Vector const &u0)
 Constructor. More...
 
void setApproximateQuadratureMatrix (SDCTimeGrid::RealMatrix const &Shat_)
 Sets the approximate quadrature matrix. More...
 
void setInitialValue (Vector const &u0, double decayRate=0.0)
 Sets \( u \gets u + e^{-(t-t_0)d} (u_0-u) \). More...
 

Constructor & Destructor Documentation

◆ Sdc()

template<class Vector >
Kaskade::Sdc< Vector >::Sdc ( SDCTimeGrid const &  ts_,
Vector const &  u0 
)
inline

Constructor.

Definition at line 1065 of file sdc.hh.

Member Function Documentation

◆ setApproximateQuadratureMatrix()

template<class Vector >
void Kaskade::Sdc< Vector >::setApproximateQuadratureMatrix ( SDCTimeGrid::RealMatrix const &  Shat_)
inline

Sets the approximate quadrature matrix.

Parameters
Shatlower triangular approximate quadrature matrix

Definition at line 1074 of file sdc.hh.

◆ setInitialValue()

template<class Vector >
void Kaskade::Sdc< Vector >::setInitialValue ( Vector const &  u0,
double  decayRate = 0.0 
)
inline

Sets \( u \gets u + e^{-(t-t_0)d} (u_0-u) \).

Definition at line 1082 of file sdc.hh.


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