Provides routines to solve systems of ODEs.
An interface for the (linear) left-hand-side of an ODE being solved by quasilinearisation.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=r8), | intent(in), | dimension(:) | :: | u | The state vector for the system of differential equations  | 
  
An interface for the (nonlinear) right-hand-side of an ODE being solved by quasilinearisation.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=r8), | intent(in), | dimension(:,:) | :: | u | The state vector for the system of differential equations, and its derivatives. Column represents the derivative.  | 
  
An interface for the product of the Jacobian of the (nonlinear) right-hand-side of an ODE and another vector.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| 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.  | 
  
An interface for a function evaluating the derivative of the state vector.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| 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  | 
  
An interface for a preconditioner to be used with the quasilinearisation ODE solver.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| 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   | 
  
|
| real(kind=r8), | intent(in), | dimension(:) | :: | rhs | The right hand side of the linear system being preconditioned.  | 
  
The result of applying the preconditioner.
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.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| 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   | 
  
||
| 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.  | 
  
||
| 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
   | 
  
|
| 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   | 
  
||
| 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.  |