module Arkode_bbd:sig
..end
Parallel band-block-diagonal preconditioners for ARKODE (requires MPI).
type[< `ARKStep | `MRIStep ]
parallel_session =(Nvector_parallel.data, Nvector_parallel.kind,
[< `ARKStep | `MRIStep ] as 'a)
session
Alias for sessions based on parallel nvectors.
type[< `ARKStep | `MRIStep ]
parallel_preconditioner =(Nvector_parallel.data, Nvector_parallel.kind,
[< `ARKStep | `MRIStep ] as 'a)
Arkode.Spils.preconditioner
Alias for preconditioners based on parallel nvectors.
typebandwidths =
ArkodeBbdTypes.bandwidths
= {
|
mudq : |
(* | Upper half-bandwidth for the difference quotient Jacobian approximation. | *) |
|
mldq : |
(* | Lower half-bandwidth for the difference quotient Jacobian approximation. | *) |
|
mukeep : |
(* | Upper half-bandwidth for the retained banded approximate Jacobian block. | *) |
|
mlkeep : |
(* | Lower half-bandwidth for the retained banded approximate Jacobian block. | *) |
}
The bandwidths for the difference quotient Jacobian operation.
typelocal_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.
typecomm_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.
val prec_left : ?dqrely:float ->
bandwidths ->
?comm:comm_fn ->
local_fn ->
[< `ARKStep | `MRIStep ] 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
.
NB: Band-Block-Diagonal preconditioners may not be used for problems involving a non-identity mass matrix.
val prec_right : ?dqrely:float ->
bandwidths ->
?comm:comm_fn ->
local_fn ->
[< `ARKStep | `MRIStep ] 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
.
val prec_both : ?dqrely:float ->
bandwidths ->
?comm:comm_fn ->
local_fn ->
[< `ARKStep | `MRIStep ] 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
.
val reinit : [< `ARKStep | `MRIStep ] 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.
val get_work_space : [< `ARKStep | `MRIStep ] parallel_session -> int * int
Returns the sizes of the real and integer workspaces used by the BBD preconditioner.
real_size
, integer_size
)val get_num_gfn_evals : [< `ARKStep | `MRIStep ] 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.