ode_solvers_mod Module

Provides routines to solve systems of ODEs.


Uses

  • module~~ode_solvers_mod~~UsesGraph module~ode_solvers_mod ode_solvers_mod penf penf module~ode_solvers_mod->penf module~nitsol_mod nitsol_mod module~ode_solvers_mod->module~nitsol_mod logger_mod logger_mod module~ode_solvers_mod->logger_mod iso_fortran_env iso_fortran_env module~ode_solvers_mod->iso_fortran_env module~nitsol_mod->iso_fortran_env

Used by

  • module~~ode_solvers_mod~~UsedByGraph module~ode_solvers_mod ode_solvers_mod module~asymmetric_plume_mod asymmetric_plume_mod module~asymmetric_plume_mod->module~ode_solvers_mod module~static_plume_mod static_plume_mod module~static_plume_mod->module~ode_solvers_mod module~plume_mod plume_mod module~plume_mod->module~ode_solvers_mod

Contents


Abstract Interfaces

abstract interface

  • public function L_intr(u)

    An interface for the (linear) left-hand-side of an ODE being solved by quasilinearisation.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in), dimension(:):: u

    The state vector for the system of differential equations

    Return Value real(kind=r8), dimension(size(u))

abstract interface

  • public function f_intr(u)

    An interface for the (nonlinear) right-hand-side of an ODE being solved by quasilinearisation.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in), dimension(:,:):: u

    The state vector for the system of differential equations, and its derivatives. Column represents the derivative.

    Return Value real(kind=r8), dimension(size(u,1))

abstract interface

  • public function jac_intr(u, du)

    An interface for the product of the Jacobian of the (nonlinear) right-hand-side of an ODE and another vector.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in), dimension(:,:):: u

    The state vector for the system of differential equations, and its derivatives, for which the Jacobian should be evaluated. Column represents the derivative.

    real(kind=r8), intent(in), dimension(:,:):: du

    The state vector for the system of differential equations, and its derivatives, which the Jacobian operates on. Column represents the derivative.

    Return Value real(kind=r8), dimension(size(u,1))

abstract interface

  • public function diff_intr(u, n)

    An interface for a function evaluating the derivative of the state vector.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in), dimension(:):: u

    The state vector for the system of differential equations

    integer, intent(in) :: n

    The order of the derivative to take

    Return Value real(kind=r8), dimension(size(u))

abstract interface

  • public function pre_intr(v, u, L, f, fcur, rhs)

    An interface for a preconditioner to be used with the quasilinearisation ODE solver.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in), dimension(:):: v

    The vector to be preconditioned.

    real(kind=r8), intent(in), dimension(:,:):: u

    The current state vector for the system of differential equations, and its derivatives. Column represents the derivative.

    procedure(L_intr) :: L

    The linear, left-hand-side of the ODE being solved.

    procedure(f_intr) :: f

    The nonlinear, right-hand-side of the ODE being solved.

    real(kind=r8), intent(in), dimension(:):: fcur

    The result of f(u)

    real(kind=r8), intent(in), dimension(:):: rhs

    The right hand side of the linear system being preconditioned.

    Return Value real(kind=r8), dimension(size(v))

    The result of applying the preconditioner.


Subroutines

public subroutine quasilinear_solve(L, f, jac_prod, solution, order, resid_norm, flag, info, tol, precond, differentiate, iter_max, gmres_iter_max, krylov_dim)

Author
Chris MacMackin
Date
March 2017

This is an iterative method to solve nonlinear systems of ODEs. Consider For a domain from , boundary conditions are specified by Here, is an th order ordinary differential operator, and are nonlinear functions of and its first derivatives.

Read more…

Arguments

Type IntentOptional AttributesName
procedure(L_intr) :: L

A function providing the linear, left-hand-side of the ODE being solved.

procedure(f_intr) :: f

A function providing the nonlinear, right-hand-side of the ODE being solved.

procedure(jac_intr) :: jac_prod

A function providing the product of the Jacobian of the nonlinear, right-hand-side of the ODE being solved and another vector.

real(kind=r8), intent(inout), dimension(:):: solution

On input, an estimate of the solution to the ODE. On output, the actual solution.

integer, intent(in) :: order

The order of the derivative taken by L

real(kind=r8), intent(out) :: resid_norm

Norm of the residual of the final solution.

integer, intent(out) :: flag

Status flag indicating whether the iterations ended succesfully.

Read more…
integer, intent(out), optional dimension(5):: info

Array containing various outputs; see above

real(kind=r8), intent(in), optional :: tol

The required reduction in the solution residual. Default is size(solution) * 1e-8.

procedure(pre_intr), optional :: precond

A right-preconditioner which may be used to improve convergence of the solution.

procedure(diff_intr), optional :: differentiate

A procedure which will evaluate the nth derivative of the state vector, when n is less than order.

integer, intent(in), optional :: iter_max

Maximum allowable number of quasilinearised iterations. Default is 15.

integer, intent(in), optional :: gmres_iter_max

Maximum allowable number of GMRES iterations. Default is 1000.

integer, intent(in), optional :: krylov_dim

Maximum Krylov subspace dimension; default 10. Larger values will allow for faster convergence (and in some cases be the difference between whether or not convergence is possible), but require more memory.