Module Sundials_LinearSolver.Iterative.Algorithms

module Algorithms: sig .. end

Low-level routines on arrays.


Scaled Preconditioned Iterative Linear Solvers routines.

Global constants and general purpose solver routines.

val qr_fact : int -> Sundials.RealArray2.t -> Sundials.RealArray.t -> bool -> unit

Performs a QR factorization of a Hessenberg matrix. The call qr_fact n h q factored, where h is the n+1 by n Hessenberg matrix (stored row-wise), q stores the computed Givens rotation, and factored=true indicates that the first n-1 columns of h have already been factored. The computed Givens rotation has the form $\begin{bmatrix} c & -s \\ s & c \end{bmatrix}$ . It is stored in the 2n elements of q as [|c; s; c; s; ...; c; s|].

val qr_sol : int ->
Sundials.RealArray2.t -> Sundials.RealArray.t -> Sundials.RealArray.t -> unit

Solve the linear least squares problem. In qr_sol n h q b, h and q are, respectively, the upper triangular factor $R$ of the original Hessenberg matrix and Q the Givens rotations used to factor it—both computed by Sundials_LinearSolver.Iterative.Algorithms.qr_fact. The function computes the n+1 elements of b to solve $Rx = Qb$.

val modified_gs : ('d, 'k) Nvector.t array -> Sundials.RealArray2.t -> int -> int -> float

Performs a modified Gram-Schmidt orthogonalization. In modified_gs v h k p,

The vector v[k] is orthogonalized against the p unit vectors at v.{k-1}, v.{k-2}, ..., v.{k-p}. The matrix h must be allocated row-wise so that the (i,j)th entry is h.{i}.{j}. The inner products are computed, $\mathtt{h.\{}i\mathtt{, k-1\}} = \mathtt{v.\{}i\mathtt{\}} \cdot \mathtt{v.\{k\}}$ , for $i=\max(0, \mathtt{k}-\mathtt{p})\ldots \mathtt{k}-1$ . The orthogonalized v.{k} is not normalized and is stored over the old v.{k}. The function returns the Euclidean norm of the orthogonalized vector.

val classical_gs : ('d, 'k) Nvector.t array ->
Sundials.RealArray2.t ->
int -> int -> Sundials.RealArray.t -> ('d, 'k) Nvector.t array -> float

Performs a classical Gram-Schmidt orthogonalization. In classical_gs v h k p s temp,

The vector v[k] is orthogonalized against the p unit vectors at v.{k-1}, v.{k-2}, ..., v.{k-p}. The matrix h must be allocated row-wise so that the (i,j)th entry is h.{i}.{j}. The inner products are computed, $\mathtt{h.\{}i\mathtt{, k-1\}} = \mathtt{v.\{}i\mathtt{\}} \cdot \mathtt{v.\{k\}}$ , for $i=\max(0, \mathtt{k}-\mathtt{p})\ldots \mathtt{k}-1$ . The orthogonalized v.{k} is not normalized and is stored over the old v.{k}. The function returns the Euclidean norm of the orthogonalized vector.