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:Sundials_LinearSolver.Iterative.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:Sundials_LinearSolver.Iterative.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 ->
    Sundials_LinearSolver.Iterative.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 ->
    Sundials_LinearSolver.Iterative.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
end