module Spils:sig..end
Scaled Preconditioned Iterative Linear Solvers.
include Sundials_LinearSolver.Iterative
include Arkode.Spils
val prec_none : ('d, 'k, Arkode.mristep) preconditionerNo preconditioning.
val prec_left : ?setup:'d prec_setup_fn ->
'd prec_solve_fn -> ('d, 'k, Arkode.mristep) preconditionerLeft preconditioning. $(P^{-1}A)x = P^{-1}b$ .
val prec_right : ?setup:'d prec_setup_fn ->
'd prec_solve_fn -> ('d, 'k, Arkode.mristep) preconditionerRight preconditioning. $(AP^{-1})Px = b$ .
val prec_both : ?setup:'d prec_setup_fn ->
'd prec_solve_fn -> ('d, 'k, Arkode.mristep) preconditionerLeft and right preconditioning. $(P_L^{-1}AP_R^{-1})P_Rx = P_L^{-1}b$
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_solverCreate an MRIStep-specific linear solver from a generic iterative linear solver.
The jac_times_rhs argument specifies an alternative right-hand-side
function for use in the internal Jacobian-vector product difference
quotient approximation. It is incorrect to specify both this argument
and jac_times_vec.
NB: a jac_times_setup_fn is not supported in
Config.sundials_version < 3.0.0.
NB: a jac_times_rhs function is not supported in
Config.sundials_version < 5.3.0.
val set_jac_eval_frequency : ('d, 'k) Arkode.MRIStep.session -> int -> unitSets the maximum number of time steps to wait before recomputation of the Jacobian or recommendation to update the preconditioner. If the integer argument is less than or equal to 0, a default value of 50 is used.
val set_linear_solution_scaling : ('d, 'k) Arkode.MRIStep.session -> bool -> unitEnables or disables scaling of the linear system solution to account for a change in $\gamma$ in the linear system. Linear solution scaling is enabled by default when a matrix-based linear solver is attached.
val set_eps_lin : ('d, 'k) Arkode.MRIStep.session -> float -> unitSets the factor by which the Krylov linear solver's convergence test constant is reduced from the Newton iteration test constant. This factor must be >= 0; passing 0 specifies the default (0.05).
val set_ls_norm_factor : ('d, 'k) Arkode.MRIStep.session -> float -> unitSets the factor for converting from the integrator tolerance (WRMS norm) to the linear solver tolerance (L2 norm). That is, $\mathit{tol}_{\mathsf{L2}} = \mathit{fact}\cdot\mathit{tol}_{\mathsf{WRMS}}$ . The given value is used directly if it is greater than zero. If it is zero (the default), then the square root of the state vector length is used. If it is less than zero, then the square root of the dot product of a state vector full of ones with itself is used.
val get_work_space : ('d, 'k) Arkode.MRIStep.session -> int * intReturns the sizes of the real and integer workspaces used by the spils linear solver.
real_size, integer_size)val get_num_lin_iters : ('d, 'k) Arkode.MRIStep.session -> intReturns the cumulative number of linear iterations.
val get_num_lin_conv_fails : ('d, 'k) Arkode.MRIStep.session -> intReturns the cumulative number of linear convergence failures.
val get_num_prec_evals : ('d, 'k) Arkode.MRIStep.session -> intReturns the cumulative number of calls to the setup function with
jok=false.
val get_num_prec_solves : ('d, 'k) Arkode.MRIStep.session -> intReturns the cumulative number of calls to the preconditioner solve function.
val get_num_jtsetup_evals : ('d, 'k) Arkode.MRIStep.session -> intReturns the cumulative number of calls to the Jacobian-vector setup function.
val get_num_jtimes_evals : ('d, 'k) Arkode.MRIStep.session -> intReturns the cumulative number of calls to the Jacobian-vector function.
val get_num_lin_rhs_evals : ('d, 'k) Arkode.MRIStep.session -> intReturns the number of calls to the right-hand side callback for finite difference Jacobian-vector product approximation. This counter is only updated if the default difference quotient function is used.
The Arkode.MRIStep.init and Arkode.MRIStep.reinit functions are the preferred way to set or
change preconditioner functions. These low-level functions are
provided for experts who want to avoid resetting internal counters
and other associated side-effects.
val set_preconditioner : ('d, 'k) Arkode.MRIStep.session ->
?setup:'d prec_setup_fn -> 'd prec_solve_fn -> unitChange the preconditioner functions.
val set_jac_times : ('d, 'k) Arkode.MRIStep.session ->
?jac_times_setup:'d jac_times_setup_fn -> 'd jac_times_vec_fn -> unitChange the Jacobian-times-vector function.
val clear_jac_times : ('d, 'k) Arkode.MRIStep.session -> unitRemove a Jacobian-times-vector function and use the default implementation.