MathTL
 All Classes Functions Variables Typedefs Enumerations
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Attributes
MathTL::WMethod< VECTOR > Class Template Reference

#include <w_method.h>

Inheritance diagram for MathTL::WMethod< VECTOR >:
MathTL::OneStepScheme< VECTOR > MathTL::ROWMethod< VECTOR >

List of all members.

Public Types

enum  Method {
  GRK4T, RODAS3, ROS2, ROS3,
  ROWDA3, ROS3P, ROS3Pw, ROSI2P2,
  ROSI2PW, RODAS, RODASP
}

Public Member Functions

 WMethod (const Method method, const WMethodStageEquationHelper< VECTOR > *stage_equation_solver)
virtual ~WMethod ()
void increment (const AbstractIVP< VECTOR > *ivp, const double t_m, const VECTOR &u_m, const double tau, VECTOR &u_mplus1, VECTOR &error_estimate, const double tolerance=1e-2) const
int order () const
void set_preprocessor (WMethodPreprocessRHSHelper< VECTOR > *pp)

Static Public Member Functions

static void transform_coefficients (const LowerTriangularMatrix< double > &Alpha, const LowerTriangularMatrix< double > &Gamma, const Vector< double > &b, const Vector< double > &bhat, LowerTriangularMatrix< double > &A, LowerTriangularMatrix< double > &C, Vector< double > &m, Vector< double > &e, Vector< double > &alpha_vector, Vector< double > &gamma_vector)
static void check_order_conditions (const LowerTriangularMatrix< double > &Alpha, const LowerTriangularMatrix< double > &Gamma, const Vector< double > &b, const Vector< double > &bhat, const bool wmethod=false)

Public Attributes

LowerTriangularMatrix< double > A
 A=(a_{i,j})
LowerTriangularMatrix< double > C
 C=(c_{i,j})
Vector< double > m
 m=(m_i)
Vector< double > e
 coefficients for the lower order error estimator e=(e_i)
Vector< double > alpha_vector
 (alpha_i)
Vector< double > gamma_vector
 (gamma_i)

Protected Attributes

const
WMethodStageEquationHelper
< VECTOR > * 
stage_equation_helper
 instance of the stage equation helper
WMethodPreprocessRHSHelper
< VECTOR > * 
preprocessor
 instance of the RHS preprocessor
int p
 consistency order

Detailed Description

template<class VECTOR>
class MathTL::WMethod< VECTOR >

The following class is an abstract base for an s-stage linearly implicit method (a W-method) for the numerical solution of (abstract) nonautonomous initial value problems of the form

u'(t) = f(t, u(t)), u(0) = u_0.

The numerical method solves the following s linear stage equations:

((*{i,i})^{-1}I - T) u_i = f(t_m + * , u^{(m)} + {j=1}^{i-1} a_{i,j} * u_j) + {j=1}^{i-1} {c_{i,j}}{} * u_j + * * g

This form arises when reformulating the original stage equations

(I-*{i,i}*T)k_i = f(t_m + * , u^{(m)} + * {j=1}^{i-1} {i,j} * k_j) + * T{j=1}^{i-1} {i,j} * k_j + * * g

by introducing the new variables u_i={j=1}^{i-1}{i,j}k_j, C=(c_{i,j})=diag(Gamma)^{-1}-Gamma^{-1}, A=(a_{i,j})=(alpha_{i,j})Gamma^{-1}. After this transformation, quite a few multiplications can be saved. The original stages and the new solution can be recovered via

k_i = 1/gamma_{i,i}*u_i - {j=1}^{i-1}c_{i,j}u_j

u^{(m+1)} = u^{(m)} + * {i=1}^s b_i*k_i = u^{(m)} + {i=1}^s m_i*u_i,

where (m_1,...,m_s)=(b_1,...,b_s)*Gamma^{-1}. The values gamma_{i,i} are stored in the diagonal of C. We assume that = {j=1}^{i-1}{i,j}, = {j=1}^i{i,j}. T is the exact Jacobian f_v(t_m,u^{(m)}) (ROW-method) or an approximation of it (W-method). g is the exact derivative f_t(t_m,u^{(m)}) (ROW-method) or an approximation of it (W-method). In a general W-method, one will often choose g=0.

The template parameter VECTOR stands for an element of the Hilbert/Banach space X under consideration. For finite-dimensional problems, this will almost always be a vector class like Vector<double>. In the infinite-dimensional case, VECTOR models elements of X as (finite) linear combinations of a wavelet basis.

All methods under consideration provide embedded lower order error estimators, of the form

uhat^{(m+1)} = u^{(m)} + * {i=1}^s bhat_i*k_i


Member Enumeration Documentation

template<class VECTOR>
enum MathTL::WMethod::Method

enum type for the different builtin methods

References: [GRK4T] Kaps, Rentrop: Generalized Runge-Kutta methods of order four with stepsize control for stiff ordinary differential equations, Numer. Math. 33(1979), 55-68 [ROS3], [RODAS3] Blom, Carmichael, Potra, Sandu, Spee, Verwer: Benchmarking Stiff ODE Solvers for Atmospheric Chemistry Problems II: Rosenbrock Solvers, Atmos. Environ. 31(1997), 3459-3472 [ROS2] Blom, Hundsdorfer, Spee, Verwer: A Second-Order Rosenbrock Method Applied to Photochemical Dispersion Problems, SIAM J. Sci. Comput. 20(1999), 1456-1480 [ROWDA3] Roche: Rosenbrock Methods for Differential Algebraic Equations, Numer. Math. 52(1988), 45-63 [ROS3P] Lang, Verwer: ROS3P - An accurate third-order Rosenbrock solver designed for parabolic problems, BIT 41(2001), 731--738 [ROS3Pw] Rang, Angermann: Creating new Rosenbrock methods with Maple, Mathematik-Bericht 2003/5, Institut fuer Mathematik, TU Clausthal [ROSI2P2], [ROSI2PW] Rang, Angermann: New Rosenbrock methods of order 3 for PDAEs of index 2 [RODAS] Hairer, Wanner: Solving Ordinairy Differential Equations II Springer Series in Computational Mathematics, Springer 1991 [RODASP] Steinebach: Order-reduction of ROW-methods for DAEs and method of lines applications Technical report 1741, TH Darmstadt, 1995


Constructor & Destructor Documentation

template<class VECTOR >
MathTL::WMethod< VECTOR >::WMethod ( const Method  method,
const WMethodStageEquationHelper< VECTOR > *  stage_equation_solver 
)

constructor from some builtin methods and a given stage equation solver

template<class VECTOR>
virtual MathTL::WMethod< VECTOR >::~WMethod ( ) [inline, virtual]

virtual destructor


Member Function Documentation

template<class VECTOR >
void MathTL::WMethod< VECTOR >::check_order_conditions ( const LowerTriangularMatrix< double > &  Alpha,
const LowerTriangularMatrix< double > &  Gamma,
const Vector< double > &  b,
const Vector< double > &  bhat,
const bool  wmethod = false 
) [static]

helper function to check the algebraic order conditions of a (RO)W-method

template<class VECTOR >
void MathTL::WMethod< VECTOR >::increment ( const AbstractIVP< VECTOR > *  ivp,
const double  t_m,
const VECTOR &  u_m,
const double  tau,
VECTOR &  u_mplus1,
VECTOR &  error_estimate,
const double  tolerance = 1e-2 
) const [virtual]

increment function u^{(m)} -> u^{(m+1)}, also returns a local error estimator

Implements MathTL::OneStepScheme< VECTOR >.

template<class VECTOR>
int MathTL::WMethod< VECTOR >::order ( ) const [inline, virtual]

consistency/convergence order

Implements MathTL::OneStepScheme< VECTOR >.

template<class VECTOR>
void MathTL::WMethod< VECTOR >::set_preprocessor ( WMethodPreprocessRHSHelper< VECTOR > *  pp) [inline]

set preprocessor for the right-hand side

template<class VECTOR >
void MathTL::WMethod< VECTOR >::transform_coefficients ( const LowerTriangularMatrix< double > &  Alpha,
const LowerTriangularMatrix< double > &  Gamma,
const Vector< double > &  b,
const Vector< double > &  bhat,
LowerTriangularMatrix< double > &  A,
LowerTriangularMatrix< double > &  C,
Vector< double > &  m,
Vector< double > &  e,
Vector< double > &  alpha_vector,
Vector< double > &  gamma_vector 
) [static]

helper function to compute the data {A, C, m, e, alpha_vector, gamma_vector} from the original data {Alpha, Gamma, b, bhat}, cf. Hairer/Wanner II.


The documentation for this class was generated from the following files:
 All Classes Functions Variables Typedefs Enumerations