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

spectral time grid for defect correction methods with Lobatto points on \( [a,b]\) More...

#include <sdc.hh>

Detailed Description

spectral time grid for defect correction methods with Lobatto points on \( [a,b]\)

For Lobatto points, \( t_0 = a \) and \( t_n = b \) holds. Since \( t_0 \) is actually used as quadrature node, there are \( n+1 \) quadrature nodes in use. Thus, Lobatto quadrature as defined here is exact for polynomials of degree up to \( 2n-1 \).

Note that collocation with Lobatto points is A-stable but not L-stable, and therefore SDC methods on Lobatto grids may be a suboptimal choice for highly stiff problems (e.g., parabolic equations or Dirichlet b.c. realized as quadratic penalty). Consider using RadauTimeGrid in these cases.

Definition at line 373 of file sdc.hh.

Inheritance diagram for Kaskade::LobattoTimeGrid:
Kaskade::SDCTimeGrid

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

 LobattoTimeGrid (int n, double a, double b)
 constructs a Lobatto grid with \( n+1 \) points on \( [a,b] \). More...
 
virtual RealVector const & points () const
 Time points in the time step. More...
 
virtual RealMatrix const & integrationMatrix () const
 Integration matrix \( S \). More...
 
virtual RealMatrix const & differentiationMatrix () const
 Differentiation matrix \( D \). More...
 
virtual void refine (RealMatrix &p)
 perform refinement of the grid, filling the prolongation matrix More...
 
virtual RealMatrix interpolate (RealVector const &x) const
 Compute interpolation coefficients. More...
 
double point (int k) const
 Time points in the time step. 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
inherited

The type used for real vectors.

Definition at line 73 of file sdc.hh.

Constructor & Destructor Documentation

◆ LobattoTimeGrid()

Kaskade::LobattoTimeGrid::LobattoTimeGrid ( int  n,
double  a,
double  b 
)

constructs a Lobatto grid with \( n+1 \) points on \( [a,b] \).

This may throw LinearAlgebraException.

Member Function Documentation

◆ differentiationMatrix()

virtual RealMatrix const & Kaskade::LobattoTimeGrid::differentiationMatrix ( ) const
inlinevirtual

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 \]

Implements Kaskade::SDCTimeGrid.

Definition at line 393 of file sdc.hh.

◆ integrationMatrix()

virtual RealMatrix const & Kaskade::LobattoTimeGrid::integrationMatrix ( ) const
inlinevirtual

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 \]

Implements Kaskade::SDCTimeGrid.

Definition at line 388 of file sdc.hh.

◆ interpolate()

virtual RealMatrix Kaskade::LobattoTimeGrid::interpolate ( RealVector const &  x) const
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. \]

Implements Kaskade::SDCTimeGrid.

◆ point()

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

Time points in the time step.

This is a convenice function.

Definition at line 94 of file sdc.hh.

◆ points()

virtual RealVector const & Kaskade::LobattoTimeGrid::points ( ) const
inlinevirtual

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.

Implements Kaskade::SDCTimeGrid.

Definition at line 383 of file sdc.hh.

◆ refine()

virtual void Kaskade::LobattoTimeGrid::refine ( RealMatrix p)
virtual

perform refinement of the grid, filling the prolongation matrix

This may throw LinearAlgebraException.

Implements Kaskade::SDCTimeGrid.


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