TimeStepping Class Reference

To solve time stepping problems, i.e. systems of linear ordinary differential equations of the form [A2]{y"} + [A1]{y'} + [A0]{y} = {b}. More...

#include <TimeStepping.h>

Public Member Functions

 TimeStepping ()
 Default constructor.
 
 TimeStepping (TimeScheme s, real_t time_step=theTimeStep, real_t final_time=theFinalTime)
 Constructor using time discretization data.
 
 ~TimeStepping ()
 Destructor.
 
void set (TimeScheme s, real_t time_step=theTimeStep, real_t final_time=theFinalTime)
 Define data of the differential equation or system.
 
void setLinearSolver (LinearSolver &ls)
 Set reference to LinearSolver instance.
 
void setPDE (Equa &eq, bool nl=false)
 Define partial differential equation to solve.
 
void setRK4RHS (Vect< real_t > &f)
 Set intermediate right-hand side vector for the Runge-Kutta method.
 
void setRK3_TVDRHS (Vect< real_t > &f)
 Set intermediate right-hand side vector for the TVD Runge-Kutta 3 method.
 
void setInitial (Vect< real_t > &u)
 Set initial condition for the system of differential equations.
 
void setInitial (Vect< real_t > &u, Vect< real_t > &v)
 Set initial condition for a system of differential equations.
 
void setInitialRHS (Vect< real_t > &f)
 Set initial RHS for a system of differential equations when the used scheme requires it.
 
void setRHS (Vect< real_t > &b)
 Set right-hand side vector.
 
void setRHS (string exp)
 Set right-hand side as defined by a regular expression.
 
void setBC (Vect< real_t > &u)
 Set vector containing boundary condition to enforce.
 
void setBC (int code, string exp)
 Set boundary condition as defined by a regular expression.
 
void setNewmarkParameters (real_t beta, real_t gamma)
 Define parameters for the Newmark scheme.
 
void setConstantMatrix ()
 Say that matrix problem is constant.
 
void setNonConstantMatrix ()
 Say that matrix problem is variable.
 
void setLinearSolver (Iteration s=DIRECT_SOLVER, Preconditioner p=DIAG_PREC)
 Set linear solver data.
 
void setNLTerm0 (Vect< real_t > &a0, Matrix< real_t > &A0)
 Set vectors defining a nonlinear first order system of ODEs.
 
void setNLTerm (Vect< real_t > &a0, Vect< real_t > &a1, Vect< real_t > &a2)
 Set vectors defining a nonlinear second order system of ODEs.
 
real_t runOneTimeStep ()
 Run one time step.
 
void run (bool opt=false)
 Run the time stepping procedure.
 
void Assembly (const Element &el, real_t *b, real_t *A0, real_t *A1, real_t *A2=nullptr)
 Assemble element arrays into global matrix and right-hand side.
 
void SAssembly (const Side &sd, real_t *b, real_t *A=nullptr)
 Assemble side arrays into global matrix and right-hand side.
 
LinearSolvergetLSolver ()
 Return LinearSolver instance.
 

Detailed Description

To solve time stepping problems, i.e. systems of linear ordinary differential equations of the form [A2]{y"} + [A1]{y'} + [A0]{y} = {b}.

Author
Rachid Touzani

Features:

  • The system may be first or second order (first and/or second order time derivatives
  • The following time integration schemes can be used:
    • For first order systems: The following schemes are implemented Forward Euler (value: FORWARD_EULER)
      Backward Euler (value: BACKWARD_EULER)
      Crank-Nicolson (value: CRANK_NICOLSON)
      Heun (value: HEUN)
      2nd Order Adams-Bashforth (value: AB2)
      4-th order Runge-Kutta (value: RK4)
      2nd order Backward Differentiation Formula (value: BDF2)
    • For second order systems: The following schemes are implemented Newmark (value: NEWMARK)

Constructor & Destructor Documentation

◆ TimeStepping()

TimeStepping ( TimeScheme  s,
real_t  time_step = theTimeStep,
real_t  final_time = theFinalTime 
)

Constructor using time discretization data.

Parameters
[in]sChoice of the scheme: To be chosen in the enumerated variable TimeScheme (see the presentation of the class)
[in]time_stepValue of the time step. This value will be modified if an adaptive method is used. The default value for this parameter if the value given by the global variable theTimeStep
[in]final_timeValue of the final time (time starts at 0). The default value for this parameter is the value given by the global variable theFinalTime

Member Function Documentation

◆ Assembly()

void Assembly ( const Element el,
real_t *  b,
real_t *  A0,
real_t *  A1,
real_t *  A2 = nullptr 
)

Assemble element arrays into global matrix and right-hand side.

This member function is to be called from finite element equation classes

Parameters
[in]elReference to Element class
[in]bPointer to element right-hand side
[in]A0Pointer to matrix of 0-th order term (involving no time derivative)
[in]A1Pointer to matrix of first order term (involving time first derivative)
[in]A2Pointer to matrix of second order term (involving time second derivative) [Default: nullptr]

◆ run()

void run ( bool  opt = false)

Run the time stepping procedure.

Parameters
[in]optFlag to say if problem matrix is constant while time stepping (true) or not (Default value is false)
Note
This argument is not used if the time stepping scheme is explicit

◆ runOneTimeStep()

real_t runOneTimeStep ( )

Run one time step.

Returns
Value of new time step if this one is updated

◆ SAssembly()

void SAssembly ( const Side sd,
real_t *  b,
real_t *  A = nullptr 
)

Assemble side arrays into global matrix and right-hand side.

This member function is to be called from finite element equation classes

Parameters
[in]sdReference to Side class
[in]bPointer to side right-hand side
[in]APointer to matrix [Default: nullptr]

◆ set()

void set ( TimeScheme  s,
real_t  time_step = theTimeStep,
real_t  final_time = theFinalTime 
)

Define data of the differential equation or system.

Parameters
[in]sChoice of the scheme: To be chosen in the enumerated variable TimeScheme (see the presentation of the class)
[in]time_stepValue of the time step. This value will be modified if an adaptive method is used. The default value for this parameter if the value given by the global variable theTimeStep
[in]final_timeValue of the final time (time starts at 0). The default value for this parameter is the value given by the global variable theFinalTime

◆ setBC()

void setBC ( int  code,
string  exp 
)

Set boundary condition as defined by a regular expression.

Parameters
[in]codeCode for which expression is assigned
[in]expRegular expression to assign as a function of x, y, z and t

◆ setConstantMatrix()

void setConstantMatrix ( )

Say that matrix problem is constant.

This is useful if the linear system is solved by a factorization method but has no effect otherwise

◆ setInitial() [1/2]

void setInitial ( Vect< real_t > &  u)

Set initial condition for the system of differential equations.

Parameters
[in]uVector containing initial condition for the unknown
Remarks
If a second-order differential equation is to be solved, use the the same function with two initial vectors (one for the unknown, the second for its time derivative)

◆ setInitial() [2/2]

void setInitial ( Vect< real_t > &  u,
Vect< real_t > &  v 
)

Set initial condition for a system of differential equations.

Parameters
[in]uVector containing initial condition for the unknown
[in]vVector containing initial condition for the time derivative of the unknown
Note
This function can be used to provide solution at previous time step if a restarting procedure is used.
This member function is to be used only in the case of a second order system

◆ setInitialRHS()

void setInitialRHS ( Vect< real_t > &  f)

Set initial RHS for a system of differential equations when the used scheme requires it.

Giving the right-hand side at initial time is somtimes required for high order methods like Runge-Kutta

Parameters
[in]fVector containing right-hand side at initial time. This vector is helpful for high order methods
Note
This function can be used to provide solution at previous time step if a restarting procedure is used.

◆ setLinearSolver() [1/2]

void setLinearSolver ( Iteration  s = DIRECT_SOLVER,
Preconditioner  p = DIAG_PREC 
)

Set linear solver data.

Parameters
[in]sSolver identification parameter. To be chosen in the enumeration variable Iteration:
DIRECT_SOLVER, CG_SOLVER, CGS_SOLVER, BICG_SOLVER, BICG_STAB_SOLVER, GMRES_SOLVER, QMR_SOLVER [Default: DIRECT_SOLVER]
[in]pPreconditioner identification parameter. To be chosen in the enumeration variable Preconditioner:
IDENT_PREC, DIAG_PREC, ILU_PREC [Default: DIAG_PREC]
Note
The argument p has no effect if the solver is DIRECT_SOLVER

◆ setLinearSolver() [2/2]

void setLinearSolver ( LinearSolver ls)

Set reference to LinearSolver instance.

Parameters
[in]lsReference to LinearSolver instance

◆ setNewmarkParameters()

void setNewmarkParameters ( real_t  beta,
real_t  gamma 
)

Define parameters for the Newmark scheme.

Parameters
[in]betaParameter beta [Default: 0.25]
[in]gammaParameter gamma [Default: 0.5]

◆ setNLTerm()

void setNLTerm ( Vect< real_t > &  a0,
Vect< real_t > &  a1,
Vect< real_t > &  a2 
)

Set vectors defining a nonlinear second order system of ODEs.

The ODE system has the form a2(u)'' + a1(u)' + a0(u) = 0

Parameters
[in]a0Reference to Vect instance defining the 0-th order term
[in]a1Reference to Vect instance defining the first order term
[in]a2Reference to Vect instance defining the second order term

◆ setNLTerm0()

void setNLTerm0 ( Vect< real_t > &  a0,
Matrix< real_t > &  A0 
)

Set vectors defining a nonlinear first order system of ODEs.

The ODE system has the form a1(u)' + a0(u) = 0

Parameters
[in]a0Reference to Vect instance defining the 0-th order term
[in]A0Reference to Matrix instance

◆ setNonConstantMatrix()

void setNonConstantMatrix ( )

Say that matrix problem is variable.

This is useful if the linear system is solved by a factorization method but has no effect otherwise

◆ setPDE()

void setPDE ( Equa eq,
bool  nl = false 
)

Define partial differential equation to solve.

The used equation class must have been constructed using the Mesh instance

Parameters
[in]eqReference to equation instance
[in]nlToggle to say if the considered equation is linear [Default: 0] or not

◆ setRHS()

void setRHS ( string  exp)

Set right-hand side as defined by a regular expression.

Parameters
[in]expRegular expression as a function of x, y, z and t

◆ setRK3_TVDRHS()

void setRK3_TVDRHS ( Vect< real_t > &  f)

Set intermediate right-hand side vector for the TVD Runge-Kutta 3 method.

Parameters
[in]fVector containing the RHS

◆ setRK4RHS()

void setRK4RHS ( Vect< real_t > &  f)

Set intermediate right-hand side vector for the Runge-Kutta method.

Parameters
[in]fVector containing the RHS