sig
type solver_result = Success | RootsFound | StopTimeReached
type step_stats = {
num_steps : int;
actual_init_step : float;
last_step : float;
current_step : float;
current_time : float;
}
type interpolant_type = Hermite | Lagrange
type linearity = Linear of bool | Nonlinear
type 'data error_weight_fun = 'data -> 'data -> unit
type ('data, 'kind) tolerance =
SStolerances of float * float
| SVtolerances of float * ('data, 'kind) Nvector.t
| WFtolerances of 'data error_weight_fun
val default_tolerances : ('data, 'kind) tolerance
type predictor_method =
TrivialPredictor
| MaximumOrderPredictor
| VariableOrderPredictor
| CutoffOrderPredictor
type 'd nonlin_system_data = {
tcur : float;
zpred : 'd;
zi : 'd;
fi : 'd;
gamma : float;
sdata : 'd;
}
type 'd rootsfn = float -> 'd -> Sundials.RealArray.t -> unit
val no_roots : int * 'd rootsfn
type adaptivity_args = {
h1 : float;
h2 : float;
h3 : float;
e1 : float;
e2 : float;
e3 : float;
q : int;
p : int;
}
type 'd adaptivity_fn = float -> 'd -> adaptivity_args -> float
type adaptivity_params = {
ks : (float * float * float) option;
method_order : bool;
}
type 'd adaptivity_method =
PIDcontroller of adaptivity_params
| PIcontroller of adaptivity_params
| Icontroller of adaptivity_params
| ExplicitGustafsson of adaptivity_params
| ImplicitGustafsson of adaptivity_params
| ImExGustafsson of adaptivity_params
| AdaptivityFn of 'd adaptivity_fn
type 'd rhsfn = float -> 'd -> 'd -> unit
type 'd stability_fn = float -> 'd -> float
type 'd resize_fn = 'd -> 'd -> unit
type 'd postprocess_step_fn = float -> 'd -> unit
type 'd triple = 'd * 'd * 'd
type ('t, 'd) jacobian_arg =
('t, 'd) Arkode_impl.jacobian_arg = {
jac_t : float;
jac_y : 'd;
jac_fy : 'd;
jac_tmp : 't;
}
type ('d, 'k) session = ('d, 'k, Arkode.mristep) Arkode_impl.session
type 'a serial_session = (Nvector_serial.data, 'a) Arkode.MRIStep.session
constraint 'a = [> Nvector_serial.kind ]
type ('data, 'kind) linear_solver
type 'a serial_linear_solver =
(Nvector_serial.data, 'a) Arkode.MRIStep.linear_solver
constraint 'a = [> Nvector_serial.kind ]
module Dls :
sig
val dense :
?context:Sundials.Context.t ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
'a Sundials.Matrix.dense ->
(Sundials.Matrix.Dense.t, 'a, [ `Dls ])
Sundials_LinearSolver.serial_t
val lapack_dense :
?context:Sundials.Context.t ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
'a Sundials.Matrix.dense ->
(Sundials.Matrix.Dense.t, 'a, [ `Dls ])
Sundials_LinearSolver.serial_t
val band :
?context:Sundials.Context.t ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
'a Sundials.Matrix.band ->
(Sundials.Matrix.Band.t, 'a, [ `Dls ]) Sundials_LinearSolver.serial_t
val lapack_band :
?context:Sundials.Context.t ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
'a Sundials.Matrix.band ->
(Sundials.Matrix.Band.t, 'a, [ `Dls ]) Sundials_LinearSolver.serial_t
module Klu :
sig
type ordering =
Sundials_LinearSolver_impl.Klu.ordering =
Amd
| ColAmd
| Natural
val make :
?context:Sundials.Context.t ->
?ordering:ordering ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
('s, 'a) Sundials.Matrix.sparse ->
('s Sundials.Matrix.Sparse.t, 'a, [ `Dls | `Klu ])
Sundials_LinearSolver.serial_t
val reinit :
('s Sundials.Matrix.Sparse.t, [> Nvector_serial.kind ] as 'a,
[> `Klu ])
Sundials_LinearSolver.serial_t ->
('s, 'a) Sundials.Matrix.sparse -> ?nnz:int -> unit -> unit
val set_ordering :
('s Sundials.Matrix.Sparse.t, [> Nvector_serial.kind ],
[> `Klu ])
Sundials_LinearSolver.serial_t -> ordering -> unit
end
val klu :
?context:Sundials.Context.t ->
?ordering:Klu.ordering ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
('s, 'a) Sundials.Matrix.sparse ->
('s Sundials.Matrix.Sparse.t, 'a, [ `Dls | `Klu ])
Sundials_LinearSolver.serial_t
module Superlumt :
sig
type ordering =
Sundials_LinearSolver_impl.Superlumt.ordering =
Natural
| MinDegreeProd
| MinDegreeSum
| ColAmd
val make :
?context:Sundials.Context.t ->
?ordering:ordering ->
nthreads:int ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
('s, 'a) Sundials.Matrix.sparse ->
('s Sundials.Matrix.Sparse.t, 'a, [ `Dls | `Slu ])
Sundials_LinearSolver.serial_t
val set_ordering :
('s Sundials.Matrix.Sparse.t, [> Nvector_serial.kind ],
[> `Slu ])
Sundials_LinearSolver.serial_t -> ordering -> unit
end
val superlumt :
?context:Sundials.Context.t ->
?ordering:Superlumt.ordering ->
nthreads:int ->
([> Nvector_serial.kind ] as 'a) Nvector.serial ->
('s, 'a) Sundials.Matrix.sparse ->
('s Sundials.Matrix.Sparse.t, 'a, [> `Dls | `Slu ])
Sundials_LinearSolver.serial_t
type 'm jac_fn =
(Sundials.RealArray.t Common.triple, Sundials.RealArray.t)
Common.jacobian_arg -> 'm -> unit
type 'm linsys_fn =
(Sundials.RealArray.t Common.triple, Sundials.RealArray.t)
Common.jacobian_arg -> 'm -> 'm option -> bool -> float -> bool
val solver :
?jac:'m jac_fn ->
?linsys:'m linsys_fn ->
('m, Sundials.RealArray.t, [> Nvector_serial.kind ] as 'a, [> `Dls ])
Sundials.LinearSolver.t -> 'a Arkode.MRIStep.serial_linear_solver
val get_work_space :
[> Nvector_serial.kind ] Arkode.MRIStep.serial_session -> int * int
val get_num_jac_evals :
[> Nvector_serial.kind ] Arkode.MRIStep.serial_session -> int
val get_num_lin_rhs_evals :
[> Nvector_serial.kind ] Arkode.MRIStep.serial_session -> int
end
module Spils :
sig
type gramschmidt_type =
Sundials_LinearSolver_impl.Iterative.gramschmidt_type =
ModifiedGS
| ClassicalGS
val spbcgs :
?context:Sundials.Context.t ->
?maxl:int ->
('d, 'k) Nvector.t ->
('m, 'd, 'k, [ `Iter | `Spbcgs ]) Sundials_LinearSolver.t
val spfgmr :
?context:Sundials.Context.t ->
?maxl:int ->
?max_restarts:int ->
?gs_type:gramschmidt_type ->
('d, 'k) Nvector.t ->
('m, 'd, 'k, [ `Iter | `Spfgmr ]) Sundials_LinearSolver.t
val spgmr :
?context:Sundials.Context.t ->
?maxl:int ->
?max_restarts:int ->
?gs_type:gramschmidt_type ->
('d, 'k) Nvector.t ->
('m, 'd, 'k, [ `Iter | `Spgmr ]) Sundials_LinearSolver.t
val sptfqmr :
?context:Sundials.Context.t ->
?maxl:int ->
('d, 'k) Nvector.t ->
('m, 'd, 'k, [ `Iter | `Sptfqmr ]) Sundials_LinearSolver.t
val pcg :
?context:Sundials.Context.t ->
?maxl:int ->
('d, 'k) Nvector.t ->
('m, 'd, 'k, [ `Iter | `Pcg ]) Sundials_LinearSolver.t
module Algorithms :
sig
val qr_fact :
int ->
Sundials.RealArray2.t -> Sundials.RealArray.t -> bool -> unit
val qr_sol :
int ->
Sundials.RealArray2.t ->
Sundials.RealArray.t -> Sundials.RealArray.t -> unit
val modified_gs :
('d, 'k) Nvector.t array ->
Sundials.RealArray2.t -> int -> int -> float
val classical_gs :
('d, 'k) Nvector.t array ->
Sundials.RealArray2.t ->
int ->
int -> Sundials.RealArray.t -> ('d, 'k) Nvector.t array -> float
end
val set_maxl :
('m, 'd, 'k, [< `Iter | `Pcg | `Spbcgs | `Sptfqmr ])
Sundials_LinearSolver.t -> int -> unit
val set_gs_type :
('m, 'd, 'k, [< `Iter | `Spfgmr | `Spgmr ]) Sundials_LinearSolver.t ->
gramschmidt_type -> unit
val set_max_restarts :
('m, 'd, 'k, [< `Iter | `Spfgmr | `Spgmr ]) Sundials_LinearSolver.t ->
int -> unit
type preconditioning_type =
Sundials_LinearSolver_impl.Iterative.preconditioning_type =
PrecNone
| PrecLeft
| PrecRight
| PrecBoth
val set_prec_type :
('m, 'd, 'k, [> `Iter ]) Sundials_LinearSolver.t ->
preconditioning_type -> unit
val set_info_file :
('m, 'd, 'k, [> `Iter ]) Sundials_LinearSolver.t ->
?print_level:bool -> Sundials.Logfile.t -> unit
val set_print_level :
('m, 'd, 'k, [> `Iter ]) Sundials_LinearSolver.t -> bool -> unit
type 'd prec_solve_arg = {
rhs : 'd;
gamma : float;
delta : float;
left : bool;
}
type 'd prec_solve_fn =
(unit, 'd) Common.jacobian_arg -> 'd prec_solve_arg -> 'd -> unit
type 'd prec_setup_fn =
(unit, 'd) Common.jacobian_arg -> bool -> float -> bool
type ('d, 'k, 's) preconditioner =
('d, 'k, 's) Arkode_impl.SpilsTypes.preconditioner
module Banded :
sig
type bandrange = { mupper : int; mlower : int; }
val prec_left :
bandrange ->
(Nvector_serial.data, [> Nvector_serial.kind ], 's)
preconditioner
val prec_right :
bandrange ->
(Nvector_serial.data, [> Nvector_serial.kind ], 's)
preconditioner
val prec_both :
bandrange ->
(Nvector_serial.data, [> Nvector_serial.kind ], 's)
preconditioner
type ('a, 'b) serial_session =
(Nvector_serial.data, 'a, 'b) Arkode_impl.session
constraint 'a = [> Nvector_serial.kind ]
constraint 'b = [< `ARKStep | `MRIStep ]
val get_work_space :
([> Nvector_serial.kind ], [< `ARKStep | `MRIStep ])
serial_session -> int * int
val get_num_rhs_evals :
([> Nvector_serial.kind ], [< `ARKStep | `MRIStep ])
serial_session -> int
end
type 'd jac_times_setup_fn = (unit, 'd) Common.jacobian_arg -> unit
type 'd jac_times_vec_fn =
('d, 'd) Common.jacobian_arg -> 'd -> 'd -> unit
val prec_none : ('d, 'k, Arkode.mristep) preconditioner
val prec_left :
?setup:'d prec_setup_fn ->
'd prec_solve_fn -> ('d, 'k, Arkode.mristep) preconditioner
val prec_right :
?setup:'d prec_setup_fn ->
'd prec_solve_fn -> ('d, 'k, Arkode.mristep) preconditioner
val prec_both :
?setup:'d prec_setup_fn ->
'd prec_solve_fn -> ('d, 'k, Arkode.mristep) preconditioner
val solver :
('m, 'd, 'k, [> `Iter ]) Sundials.LinearSolver.t ->
?jac_times_vec:'d jac_times_setup_fn option * 'd jac_times_vec_fn ->
?jac_times_rhs:'d rhsfn ->
('d, 'k, Arkode.mristep) preconditioner ->
('d, 'k) Arkode.MRIStep.linear_solver
val set_jac_eval_frequency :
('d, 'k) Arkode.MRIStep.session -> int -> unit
val set_linear_solution_scaling :
('d, 'k) Arkode.MRIStep.session -> bool -> unit
val set_eps_lin : ('d, 'k) Arkode.MRIStep.session -> float -> unit
val set_ls_norm_factor :
('d, 'k) Arkode.MRIStep.session -> float -> unit
val get_work_space : ('d, 'k) Arkode.MRIStep.session -> int * int
val get_num_lin_iters : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_lin_conv_fails : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_prec_evals : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_prec_solves : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_jtsetup_evals : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_jtimes_evals : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_lin_rhs_evals : ('d, 'k) Arkode.MRIStep.session -> int
val set_preconditioner :
('d, 'k) Arkode.MRIStep.session ->
?setup:'d prec_setup_fn -> 'd prec_solve_fn -> unit
val set_jac_times :
('d, 'k) Arkode.MRIStep.session ->
?jac_times_setup:'d jac_times_setup_fn -> 'd jac_times_vec_fn -> unit
val clear_jac_times : ('d, 'k) Arkode.MRIStep.session -> unit
end
val matrix_embedded_solver :
(unit, 'data, 'kind, [> `MatE ]) Sundials.LinearSolver.t ->
('data, 'kind) Arkode.MRIStep.linear_solver
module InnerStepper :
sig
type ('d, 'k) t = ('d, 'k) Arkode_impl.inner_stepper
val from_arkstep :
('d, 'k) Arkode.ARKStep.session ->
('d, 'k) Arkode.MRIStep.InnerStepper.t
type 'd evolvefn = float -> float -> 'd -> unit
type fullrhs_mode = Arkode_impl.fullrhs_mode = Start | End | Other
type 'd full_rhsfn =
float ->
'd -> 'd -> Arkode.MRIStep.InnerStepper.fullrhs_mode -> unit
type 'd resetfn = float -> 'd -> unit
val make :
?context:Sundials.Context.t ->
evolve_fn:'d Arkode.MRIStep.InnerStepper.evolvefn ->
full_rhs_fn:'d Arkode.MRIStep.InnerStepper.full_rhsfn ->
?reset_fn:'d Arkode.MRIStep.InnerStepper.resetfn ->
unit -> ('d, 'k) Arkode.MRIStep.InnerStepper.t
val add_forcing :
('d, 'k) Arkode.MRIStep.InnerStepper.t ->
float -> ('d, 'k) Nvector.t -> unit
type 'd forcing_data = {
tshift : float;
tscale : float;
forcing : 'd array;
}
val get_forcing_data :
('d, 'k) Arkode.MRIStep.InnerStepper.t ->
'd Arkode.MRIStep.InnerStepper.forcing_data
end
module Coupling :
sig
type t
val nmat : Arkode.MRIStep.Coupling.t -> int
val stages : Arkode.MRIStep.Coupling.t -> int
val method_order : Arkode.MRIStep.Coupling.t -> int
val embedding_order : Arkode.MRIStep.Coupling.t -> int
val explicit_coupling_matrices :
Arkode.MRIStep.Coupling.t -> Sundials.RealArray.t array array option
val implicit_coupling_matrices :
Arkode.MRIStep.Coupling.t -> Sundials.RealArray.t array array option
val abscissae : Arkode.MRIStep.Coupling.t -> Sundials.RealArray.t
val make :
method_order:int ->
embedding_order:int ->
?explicit:Sundials.RealArray.t array array ->
?implicit:Sundials.RealArray.t array array ->
Sundials.RealArray.t -> Arkode.MRIStep.Coupling.t
type coupling_table =
KW3
| GARK_ERK33a
| GARK_ERK45a
| GARK_IRK21a
| GARK_ESDIRK34a
| GARK_ESDIRK46a
| IMEX_GARK3a
| IMEX_GARK3b
| IMEX_GARK4
val load_table :
Arkode.MRIStep.Coupling.coupling_table -> Arkode.MRIStep.Coupling.t
val mis_to_mri :
method_order:int ->
embedding_order:int ->
Arkode.ButcherTable.t -> Arkode.MRIStep.Coupling.t
val copy : Arkode.MRIStep.Coupling.t -> Arkode.MRIStep.Coupling.t
val space : Arkode.MRIStep.Coupling.t -> int * int
val write :
?logfile:Sundials.Logfile.t -> Arkode.MRIStep.Coupling.t -> unit
end
type ('d, 'k) problem
val implicit :
?nlsolver:('data, 'kind, ('data, 'kind) Arkode.MRIStep.session,
[ `Nvec ])
Sundials_NonlinearSolver.t ->
?nlsrhsfn:'data rhsfn ->
?lsolver:('data, 'kind) Arkode.MRIStep.linear_solver ->
?linearity:linearity ->
'data rhsfn -> ('data, 'kind) Arkode.MRIStep.problem
val explicit : 'data rhsfn -> ('data, 'kind) Arkode.MRIStep.problem
val imex :
?nlsolver:('data, 'kind, ('data, 'kind) Arkode.MRIStep.session,
[ `Nvec ])
Sundials_NonlinearSolver.t ->
?nlsrhsfn:'data rhsfn ->
?lsolver:('data, 'kind) Arkode.MRIStep.linear_solver ->
?linearity:linearity ->
fsi:'data rhsfn ->
fse:'data rhsfn -> unit -> ('data, 'kind) Arkode.MRIStep.problem
val init :
?context:Sundials.Context.t ->
('data, 'kind) Arkode.MRIStep.problem ->
('data, 'kind) tolerance ->
('data, 'kind) Arkode.MRIStep.InnerStepper.t ->
?coupling:Arkode.MRIStep.Coupling.t ->
slowstep:float ->
?roots:int * 'data rootsfn ->
float ->
('data, 'kind) Nvector.t -> ('data, 'kind) Arkode.MRIStep.session
val evolve_normal :
('d, 'k) Arkode.MRIStep.session ->
float -> ('d, 'k) Nvector.t -> float * solver_result
val evolve_one_step :
('d, 'k) Arkode.MRIStep.session ->
float -> ('d, 'k) Nvector.t -> float * solver_result
val get_dky :
('d, 'k) Arkode.MRIStep.session ->
('d, 'k) Nvector.t -> float -> int -> unit
val reinit :
('d, 'k) Arkode.MRIStep.session ->
?problem:('d, 'k) Arkode.MRIStep.problem ->
?roots:int * 'd rootsfn -> float -> ('d, 'k) Nvector.t -> unit
val reset : ('d, 'k) Arkode.MRIStep.session -> float -> ('d, 'k) Nvector.t
val resize :
('d, 'k) Arkode.MRIStep.session ->
?resize_nvec:'d resize_fn -> ('d, 'k) Nvector.t -> float -> unit
val set_defaults : ('d, 'k) Arkode.MRIStep.session -> unit
val set_interpolant_type :
('d, 'k) Arkode.MRIStep.session -> interpolant_type -> unit
val set_interpolant_degree : ('d, 'k) Arkode.MRIStep.session -> int -> unit
val set_diagnostics :
?logfile:Sundials.Logfile.t -> ('d, 'k) Arkode.MRIStep.session -> unit
val clear_diagnostics : ('d, 'k) Arkode.MRIStep.session -> unit
val set_error_file :
('d, 'k) Arkode.MRIStep.session -> Sundials.Logfile.t -> unit
val set_err_handler_fn :
('d, 'k) Arkode.MRIStep.session ->
(Sundials.Util.error_details -> unit) -> unit
val clear_err_handler_fn : ('d, 'k) Arkode.MRIStep.session -> unit
val set_fixed_step : ('d, 'k) Arkode.MRIStep.session -> float -> unit
val set_max_hnil_warns : ('d, 'k) Arkode.MRIStep.session -> int -> unit
val set_max_num_steps : ('d, 'k) Arkode.MRIStep.session -> int -> unit
val set_stop_time : ('d, 'k) Arkode.MRIStep.session -> float -> unit
type 'd pre_inner_fn = float -> 'd array -> unit
val set_pre_inner_fn :
('d, 'k) Arkode.MRIStep.session -> 'd Arkode.MRIStep.pre_inner_fn -> unit
val clear_pre_inner_fn : ('d, 'k) Arkode.MRIStep.session -> unit
type 'd post_inner_fn = float -> 'd -> unit
val set_post_inner_fn :
('d, 'k) Arkode.MRIStep.session ->
'd Arkode.MRIStep.post_inner_fn -> unit
val clear_post_inner_fn : ('d, 'k) Arkode.MRIStep.session -> unit
val set_postprocess_step_fn :
('d, 'k) Arkode.MRIStep.session -> 'd postprocess_step_fn -> unit
val clear_postprocess_step_fn : ('d, 'k) Arkode.MRIStep.session -> unit
type 'd stage_predict_fn = float -> 'd -> unit
val set_stage_predict_fn :
('d, 'k) Arkode.MRIStep.session ->
'd Arkode.MRIStep.stage_predict_fn -> unit
val clear_stage_predict_fn : ('d, 'k) Arkode.MRIStep.session -> unit
val set_nonlin_conv_coef : ('d, 'k) Arkode.MRIStep.session -> float -> unit
val set_nonlin_crdown : ('d, 'k) Arkode.MRIStep.session -> float -> unit
val set_nonlin_rdiv : ('d, 'k) Arkode.MRIStep.session -> float -> unit
val set_delta_gamma_max : ('d, 'k) Arkode.MRIStep.session -> float -> unit
val set_lsetup_frequency : ('d, 'k) Arkode.MRIStep.session -> int -> unit
val set_linear : ('d, 'k) Arkode.MRIStep.session -> bool -> unit
val set_nonlinear : ('d, 'k) Arkode.MRIStep.session -> unit
val set_predictor_method :
('d, 'k) Arkode.MRIStep.session -> predictor_method -> unit
val set_max_nonlin_iters : ('d, 'k) Arkode.MRIStep.session -> int -> unit
val get_work_space : ('d, 'k) Arkode.MRIStep.session -> int * int
val get_num_steps : ('d, 'k) Arkode.MRIStep.session -> int
val get_last_step : ('d, 'k) Arkode.MRIStep.session -> float
val get_num_rhs_evals : ('d, 'k) Arkode.MRIStep.session -> int * int
val get_current_time : ('d, 'k) Arkode.MRIStep.session -> float
val get_current_state : ('d, 'k) Arkode.MRIStep.session -> 'd
val get_current_coupling :
('d, 'k) Arkode.MRIStep.session -> Arkode.MRIStep.Coupling.t
val write_coupling :
?logfile:Sundials.Logfile.t -> ('d, 'k) Arkode.MRIStep.session -> unit
val get_nonlin_system_data :
('d, 'k) Arkode.MRIStep.session -> 'd nonlin_system_data
val compute_state :
('d, 'k) Arkode.MRIStep.session ->
('d, 'k) Nvector.t -> ('d, 'k) Nvector.t -> unit
val get_current_gamma : ('d, 'k) Arkode.MRIStep.session -> float
val get_tol_scale_factor : ('d, 'k) Arkode.MRIStep.session -> float
val get_err_weights :
('d, 'k) Arkode.MRIStep.session -> ('d, 'k) Nvector.t -> unit
val get_num_lin_solv_setups : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_nonlin_solv_iters : ('d, 'k) Arkode.MRIStep.session -> int
val get_num_nonlin_solv_conv_fails : ('d, 'k) Arkode.MRIStep.session -> int
val get_nonlin_solv_stats : ('d, 'k) Arkode.MRIStep.session -> int * int
val write_session :
?logfile:Sundials.Logfile.t -> ('d, 'k) Arkode.MRIStep.session -> unit
val set_root_direction :
('d, 'k) Arkode.MRIStep.session -> Sundials.RootDirs.d array -> unit
val set_all_root_directions :
('d, 'k) Arkode.MRIStep.session -> Sundials.RootDirs.d -> unit
val set_no_inactive_root_warn : ('d, 'k) Arkode.MRIStep.session -> unit
val get_num_roots : ('d, 'k) Arkode.MRIStep.session -> int
val get_root_info :
('d, 'k) Arkode.MRIStep.session -> Sundials.Roots.t -> unit
val get_num_g_evals : ('d, 'k) Arkode.MRIStep.session -> int
val write_parameters :
?logfile:Sundials.Logfile.t -> ('d, 'k) Arkode.MRIStep.session -> unit
end