KASKADE 7 development version
Public Types | Public Member Functions | List of all members
Kaskade::SDCTimeGrid Class Referenceabstract

Abstract base class of time grids for (block) spectral defect correction methods. More...

#include <sdc.hh>

Detailed Description

Abstract base class of time grids for (block) spectral defect correction methods.

This class represents a time grid on \( [a,b]\) with \( n \) subintervals and \( n+1 \) time grid points (including \( t_0=a\), but \( t_n \le b\)). Polynomial interpolation, quadrature, and differentiation are provided using (a subset of) the grid points as nodes.

Conceptually, the collocation polynomials are defined by Lagrange (or Hermite) interpolation of their values at the initial point \( t_0 \) and \( n \) collocation points (nodes) \( t_1,\dots,t_n \). Grid points may coincide, in which case Hermite interpolation is used. Collocation conditions for \( \dot u = f(u) \) are to be satisfied only at the nodes, i.e. \( \dot u(t_i) = f(u(t_i)), \quad i=1,\dots,n\).

Definition at line 67 of file sdc.hh.

Inheritance diagram for Kaskade::SDCTimeGrid:
Kaskade::LobattoTimeGrid Kaskade::SDCGridImplementationBase Kaskade::GaussTimeGrid Kaskade::RadauTimeGrid

Public Types

typedef Dune::DynamicVector< double > RealVector
 The type used for real vectors. More...
 
typedef DynamicMatrix< double > RealMatrix
 The type used for real (dense) matrices. More...
 

Public Member Functions

virtual RealVector const & points () const =0
 Time points in the time step. More...
 
double point (int k) const
 Time points in the time step. More...
 
virtual RealMatrix const & integrationMatrix () const =0
 Integration matrix \( S \). More...
 
virtual RealMatrix const & differentiationMatrix () const =0
 Differentiation matrix \( D \). More...
 
virtual void refine (RealMatrix &p)=0
 Perform refinement of the grid, filling the prolongation matrix. More...
 
virtual RealMatrix interpolate (RealVector const &x) const =0
 Compute interpolation coefficients. More...
 

Member Typedef Documentation

◆ RealMatrix

The type used for real (dense) matrices.

Definition at line 78 of file sdc.hh.

◆ RealVector

typedef Dune::DynamicVector<double> Kaskade::SDCTimeGrid::RealVector

The type used for real vectors.

Definition at line 73 of file sdc.hh.

Member Function Documentation

◆ differentiationMatrix()

virtual RealMatrix const & Kaskade::SDCTimeGrid::differentiationMatrix ( ) const
pure virtual

Differentiation matrix \( D \).

This computes a matrix such that polynomials given by values at \( t_0, \dots,t_n\) can easily be differentiated.

The Lagrangian interpolation functions \( L_k \) are defined by \( L_k(t_i) = \delta_{ik} , \quad i=0,\dots,n\). The matrix \( D \in \mathbb{R}^{n+1\times n+1}\) contains the values

\[ D_{ik} = \dot L_k(t_i) \]

This way, if \( u \) is defined in terms of its function values \( v_i = u(t_i) \), its derivatives can be evaluated by a matrix-vector multiplication:

\[ \dot u(\tau_i) = (Dv)_i \]

Implemented in Kaskade::SDCGridImplementationBase, and Kaskade::LobattoTimeGrid.

Referenced by Kaskade::waveformRelaxationStep2().

◆ integrationMatrix()

virtual RealMatrix const & Kaskade::SDCTimeGrid::integrationMatrix ( ) const
pure virtual

Integration matrix \( S \).

This computes a matrix such that functions given by values at \( t_0,\dots,t_n\) can be easily integrated.

Interpolation is based on the nodes \(t_0\dots,t_n\). Depending on the actual implementation, the initial point \( t_0 \) might be ignored (i.e. have a quadrature weight of 0).

The Lagrangian interpolation functions \( L_k \) are defined by \( L_k(t_i) = \delta_{ik},\quad i=1,\dots,n \). The matrix \( S \in \mathbb{R}^{n\times n+1}\) contains the values

\[ S_{ik} = \int_{\tau=t_i}^{t_{i+1}} L_k(\tau) \, d\tau \]

(the leading column being zero). This way, if \( u \) is defined in terms of its function values \( v_i = u(t_i),\quad i=1,\dots,n \), the integrals can be evaluated by a matrix-vector multiplication:

\[ \int_{\tau=t_i}^{t_{i+1}} u(\tau) \, d\tau = (Sv)_i \]

Implemented in Kaskade::SDCGridImplementationBase, and Kaskade::LobattoTimeGrid.

Referenced by Kaskade::sdcExEulerIterationStep(), and Kaskade::sdcIterationStep2().

◆ interpolate()

virtual RealMatrix Kaskade::SDCTimeGrid::interpolate ( RealVector const &  x) const
pure virtual

Compute interpolation coefficients.

Returns a matrix \( w \in \mathbb{R}^{m+1\times n+1} \), such that the interpolation polynomial \( p \) to the values \( y_i \) at grid points \( t_i \) can be evaluated as

\[ p(x_i) = \sum_{j=0}^n w_{ij} y_j, \quad i=0,\dots,m. \]

Implemented in Kaskade::LobattoTimeGrid, Kaskade::RadauTimeGrid, and Kaskade::GaussTimeGrid.

◆ point()

double Kaskade::SDCTimeGrid::point ( int  k) const
inline

Time points in the time step.

This is a convenice function.

Definition at line 94 of file sdc.hh.

◆ points()

virtual RealVector const & Kaskade::SDCTimeGrid::points ( ) const
pure virtual

Time points in the time step.

The time step \( [t_0, t_n] \) contains \( n+1 \) time points \( t_i \), including the end points. Those are provided here. The time points are stored in increasing order.

Implemented in Kaskade::SDCGridImplementationBase, and Kaskade::LobattoTimeGrid.

Referenced by Kaskade::EulerSDC< Vector, TimeGrid >::integrate(), Kaskade::EulerSDC< Vector, TimeGrid >::integrateTOL(), point(), Kaskade::sdcExEulerIterationStep(), Kaskade::sdcIterationStep2(), Kaskade::Sdc< Vector >::setInitialValue(), and Kaskade::waveformRelaxationStep2().

◆ refine()

virtual void Kaskade::SDCTimeGrid::refine ( RealMatrix p)
pure virtual

Perform refinement of the grid, filling the prolongation matrix.

If the function representation is not sufficiently accurate, a finer grid of time points can be tried. This method refines the grid to \( m+1 \) time points \( s_i \), \( m>n \), and fills a prolongation matrix \( P\in \mathbb{R}^{m+1\times n+1} \) such that with \( v_i = u(t_i) \) and \( w=Pv \) it holds that \( w_i = u(s_i) \).

Different derived classes may implement this in different ways, or provide their own, more flexible ways of refining the grid.

Parameters
[out]pthe prolongation matrix

Implemented in Kaskade::LobattoTimeGrid, Kaskade::RadauTimeGrid, and Kaskade::GaussTimeGrid.


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