KASKADE 7 development version
|
spectral time grid for defect correction methods with Gauss points More...
#include <sdc.hh>
spectral time grid for defect correction methods with Gauss points
For the Gauss time grid, \( t_0 = a \) and \( t_n < b \) hold. Since \( t_0 \) is not used for quadrature, there are \( n \) actual quadrature nodes in use. Thus, Gauss quadrature as defined here is exact for polynomials of degree up to \( 2n-1 \).
Note that collocation with Gauss points is A-stable, and therefore SDC methods on Gauss grids are best for oscillatory problems (e.g., hyperbolic problems). Consider using RadauTimeGrid in other cases.
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 | |
GaussTimeGrid (int n, double a, double b) | |
constructs a Gauss grid with \( n+1 \) points on \( [a,b] \) 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... | |
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... | |
double | point (int k) const |
Time points in the time step. More... | |
Protected Attributes | |
RealVector | pts |
RealMatrix | integ |
RealMatrix | diff |
|
inherited |
|
inherited |
Kaskade::GaussTimeGrid::GaussTimeGrid | ( | int | n, |
double | a, | ||
double | b | ||
) |
constructs a Gauss grid with \( n+1 \) points on \( [a,b] \)
This may throw LinearAlgebraException.
|
inlinevirtualinherited |
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.
|
inlinevirtualinherited |
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.
|
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.
|
inlineinherited |
|
inlinevirtualinherited |
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.
|
virtual |
perform refinement of the grid, filling the prolongation matrix
This may throw LinearAlgebraException.
Implements Kaskade::SDCTimeGrid.
|
protectedinherited |
Definition at line 192 of file sdc.hh.
Referenced by Kaskade::SDCGridImplementationBase::differentiationMatrix().
|
protectedinherited |
Definition at line 191 of file sdc.hh.
Referenced by Kaskade::SDCGridImplementationBase::integrationMatrix().
|
protectedinherited |
Definition at line 190 of file sdc.hh.
Referenced by Kaskade::SDCGridImplementationBase::points().