Module Arkode_bbd

module Arkode_bbd: sig .. end
Parallel band-block-diagonal preconditioners for ARKODE (requires MPI).
Author(s): Timothy Bourke (Inria/ENS), Jun Inoue (Inria/ENS), Marc Pouzet (UPMC/ENS/Inria)
Version: 2.7.0
See sundials: Parallel band-block-diagonal preconditioner module

type parallel_session = (Nvector_parallel.data, Nvector_parallel.kind) Arkode.session 
Alias for sessions based on parallel nvectors.
type parallel_preconditioner = (Nvector_parallel.data, Nvector_parallel.kind) Arkode.Spils.preconditioner 
Alias for preconditioners based on parallel nvectors.
type bandwidths = ArkodeBbdTypes.bandwidths = {
   mudq : int; (*
Upper half-bandwidth for the difference quotient Jacobian approximation.
*)
   mldq : int; (*
Lower half-bandwidth for the difference quotient Jacobian approximation.
*)
   mukeep : int; (*
Upper half-bandwidth for the retained banded approximate Jacobian block.
*)
   mlkeep : int; (*
Lower half-bandwidth for the retained banded approximate Jacobian block.
*)
}
The bandwidths for the difference quotient Jacobian operation.
type local_fn = float -> Nvector_parallel.data -> Nvector_parallel.data -> unit 
Approximates the right-hand side function using local computations only. That is, it computes the approximation $g(t,y) \approx f_I(t,y)$. In the call gloc t y g, t is the independent variable (time), y is the input vector, and g stores the computed derivatives.

Raising Sundials.RecoverableFailure signals a recoverable error. Other exceptions signal unrecoverable errors.

See sundials: ARKLocalFn
type comm_fn = float -> Nvector_parallel.data -> unit 
Functions that perform the interprocess communication necessary for the execution of Arkode_bbd.local_fn. In the call cfn t y, t is the independent variable (time) and y is the input vector.

Raising Sundials.RecoverableFailure signals a recoverable error. Other exceptions signal unrecoverable errors.

See sundials: ARKCommFn
val prec_left : ?dqrely:float ->
bandwidths ->
?comm_fn:comm_fn ->
local_fn -> parallel_preconditioner
Left preconditioning using the Parallel Band-Block-Diagonal module. The difference quotient operation is controlled by ?dqrely, which specifies the relative increment in components of y, and Arkode_bbd.bandwidths.
See sundials: ARKBBDPrecInit
val prec_right : ?dqrely:float ->
bandwidths ->
?comm_fn:comm_fn ->
local_fn -> parallel_preconditioner
Right preconditioning using the Parallel Band-Block-Diagonal module. The difference quotient operation is controlled by ?dqrely, which specifies the relative increment in components of y, and Arkode_bbd.bandwidths.
See sundials: ARKBBDPrecInit
val prec_both : ?dqrely:float ->
bandwidths ->
?comm_fn:comm_fn ->
local_fn -> parallel_preconditioner
Preconditioning from both sides using the Parallel Band-Block-Diagonal module. The difference quotient operation is controlled by ?dqrely, the relative increment in components of y, and Arkode_bbd.bandwidths.
See sundials: ARKBBDPrecInit
val reinit : parallel_session -> ?dqrely:float -> int -> int -> unit
Reinitializes some BBD preconditioner parameters. In the call, reinit s ~dqrely:dqrely mudq mldq, dqrely is the relative increment in the components of y, and mudq and mldq are, respectively, the upper-half and lower-half bandwidths of the difference quotient Jacobian approximation.
See sundials: ARKBBDPrecReInit
val get_work_space : parallel_session -> int * int
Returns the sizes of the real and integer workspaces used by the BBD preconditioner.
Returns (real_size, integer_size)
See sundials: ARKBBDPrecGetWorkSpace
val get_num_gfn_evals : parallel_session -> int
Returns the number of calls to the right-hand side function due to finite difference banded Jacobian approximation in the setup function.
See sundials: ARKBBDPrecGetNumGfnEvals