To solve time stepping problems, i.e. systems of linear ordinary differential equations of the form [A2]{y"} + [A1]{y'} + [A0]{y} = {b}. More...
Public Member Functions | |
TimeStepping () | |
Default constructor. | |
TimeStepping (TimeScheme s, real_t time_step=theTimeStep, real_t final_time=theFinalTime) | |
Constructor using time discretization data. More... | |
~TimeStepping () | |
Destructor. | |
void | set (TimeScheme s, real_t time_step=theTimeStep, real_t final_time=theFinalTime) |
Define data of the differential equation or system. More... | |
void | setLinearSolver (LinearSolver &ls) |
Set reference to LinearSolver instance. More... | |
void | setPDE (Equa &eq, bool nl=false) |
Define partial differential equation to solve. More... | |
void | setRK4RHS (Vect< real_t > &f) |
Set intermediate right-hand side vector for the Runge-Kutta method. More... | |
void | setRK3_TVDRHS (Vect< real_t > &f) |
Set intermediate right-hand side vector for the TVD Runge-Kutta 3 method. More... | |
void | setInitial (Vect< real_t > &u) |
Set initial condition for the system of differential equations. More... | |
void | setInitial (Vect< real_t > &u, Vect< real_t > &v) |
Set initial condition for a system of differential equations. More... | |
void | setInitialRHS (Vect< real_t > &f) |
Set initial RHS for a system of differential equations when the used scheme requires it. More... | |
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. More... | |
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. More... | |
void | setNewmarkParameters (real_t beta, real_t gamma) |
Define parameters for the Newmark scheme. More... | |
void | setConstantMatrix () |
Say that matrix problem is constant. More... | |
void | setNonConstantMatrix () |
Say that matrix problem is variable. More... | |
void | setLinearSolver (Iteration s=DIRECT_SOLVER, Preconditioner p=DIAG_PREC) |
Set linear solver data. More... | |
void | setNLTerm0 (Vect< real_t > &a0, Matrix< real_t > &A0) |
Set vectors defining a nonlinear first order system of ODEs. More... | |
void | setNLTerm (Vect< real_t > &a0, Vect< real_t > &a1, Vect< real_t > &a2) |
Set vectors defining a nonlinear second order system of ODEs. More... | |
real_t | runOneTimeStep () |
Run one time step. More... | |
void | run (bool opt=false) |
Run the time stepping procedure. More... | |
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. More... | |
void | SAssembly (const Side &sd, real_t *b, real_t *A=nullptr) |
Assemble side arrays into global matrix and right-hand side. More... | |
LinearSolver & | getLSolver () |
Return LinearSolver instance. | |
To solve time stepping problems, i.e. systems of linear ordinary differential equations of the form [A2]{y"} + [A1]{y'} + [A0]{y} = {b}.
Features:
TimeStepping | ( | TimeScheme | s, |
real_t | time_step = theTimeStep , |
||
real_t | final_time = theFinalTime |
||
) |
Constructor using time discretization data.
[in] | s | Choice of the scheme: To be chosen in the enumerated variable TimeScheme (see the presentation of the class) |
[in] | time_step | Value 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_time | Value of the final time (time starts at 0). The default value for this parameter is the value given by the global variable theFinalTime |
Assemble element arrays into global matrix and right-hand side.
This member function is to be called from finite element equation classes
[in] | el | Reference to Element class |
[in] | b | Pointer to element right-hand side |
[in] | A0 | Pointer to matrix of 0-th order term (involving no time derivative) |
[in] | A1 | Pointer to matrix of first order term (involving time first derivative) |
[in] | A2 | Pointer to matrix of second order term (involving time second derivative) [Default: nullptr ] |
void run | ( | bool | opt = false | ) |
Run the time stepping procedure.
[in] | opt | Flag to say if problem matrix is constant while time stepping (true) or not (Default value is false) |
real_t runOneTimeStep | ( | ) |
Run one time step.
Assemble side arrays into global matrix and right-hand side.
This member function is to be called from finite element equation classes
[in] | sd | Reference to Side class |
[in] | b | Pointer to side right-hand side |
[in] | A | Pointer to matrix [Default: nullptr ] |
void set | ( | TimeScheme | s, |
real_t | time_step = theTimeStep , |
||
real_t | final_time = theFinalTime |
||
) |
Define data of the differential equation or system.
[in] | s | Choice of the scheme: To be chosen in the enumerated variable TimeScheme (see the presentation of the class) |
[in] | time_step | Value 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_time | Value of the final time (time starts at 0). The default value for this parameter is the value given by the global variable theFinalTime |
void setBC | ( | int | code, |
string | exp | ||
) |
Set boundary condition as defined by a regular expression.
[in] | code | Code for which expression is assigned |
[in] | exp | Regular expression to assign as a function of x , y , z and t |
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
Set initial condition for the system of differential equations.
[in] | u | Vector containing initial condition for the unknown |
Set initial condition for a system of differential equations.
[in] | u | Vector containing initial condition for the unknown |
[in] | v | Vector containing initial condition for the time derivative of the unknown |
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
[in] | f | Vector containing right-hand side at initial time. This vector is helpful for high order methods |
void setLinearSolver | ( | LinearSolver & | ls | ) |
Set reference to LinearSolver instance.
[in] | ls | Reference to LinearSolver instance |
void setLinearSolver | ( | Iteration | s = DIRECT_SOLVER , |
Preconditioner | p = DIAG_PREC |
||
) |
Set linear solver data.
[in] | s | Solver 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] | p | Preconditioner identification parameter. To be chosen in the enumeration variable Preconditioner: IDENT_PREC, DIAG_PREC, ILU_PREC [Default: DIAG_PREC ] |
Define parameters for the Newmark scheme.
[in] | beta | Parameter beta [Default: 0.25 ] |
[in] | gamma | Parameter gamma [Default: 0.5 ] |
Set vectors defining a nonlinear second order system of ODEs.
The ODE system has the form a2(u)'' + a1(u)' + a0(u) = 0
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
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
[in] | eq | Reference to equation instance |
[in] | nl | Toggle to say if the considered equation is linear [Default: 0 ] or not |
void setRHS | ( | string | exp | ) |
Set right-hand side as defined by a regular expression.
[in] | exp | Regular expression as a function of x , y , z and t |
Set intermediate right-hand side vector for the TVD Runge-Kutta 3 method.
[in] | f | Vector containing the RHS |