Runge-Kutta TVD ODE solver class.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | neq |
number of equations |
|||
integer, | public | :: | order |
order of the method |
|||
integer, | public | :: | fevals | = | 0 |
number of function evaluations |
|
integer, | public | :: | istate | = | 0 |
flag indicating the state of the integration: 1 first call for a problem, 2 subsequent call for a problem. |
|
character(len=:), | public, | allocatable | :: | msg |
error message |
Initialize rktvd
object.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(integrand) | :: | fu |
subroutine with the derivative |
|||
integer, | intent(in) | :: | neq |
number of equations |
||
integer, | intent(in) | :: | order |
order of the method (1, 2 or 3) |
Error method.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(tvdode), | intent(inout) | :: | self |
object |
||
character(len=*), | intent(in) | :: | msg |
message |
This subroutine implements the optimal 1st, 2nd and 3rd order TVD RK methods described in ICASE 97-65 (Shu, 1997). The routine was built to work similarly to LSODE.
Note
There are also 4th and 5th order methods, but they have lower CFL coefficients and are more difficult to implement. See Equation 4.15, page 44.
Todo
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rktvd), | intent(inout) | :: | self |
object |
||
real(kind=rk), | intent(inout) | :: | u(:) |
vector(neq) with the variables to integrate |
||
real(kind=rk), | intent(inout) | :: | t |
time; on return it will be the current value of (close to tout) |
||
real(kind=rk), | intent(in) | :: | tout |
time where next output is desired |
||
real(kind=rk), | intent(in) | :: | dt |
time step |
||
integer, | intent(in), | optional | :: | itask |
flag indicating the task to be performed:
1 normal integration until tout;
2 single |
type, extends(tvdode) :: rktvd !! Runge-Kutta TVD ODE solver class. contains procedure, pass(self) :: integrate => rktvd_integrate end type rktvd