Module Sundials_Matrix.ArrayBand

module ArrayBand: sig .. end

General-purpose band matrix operations on arrays.


type smu = int 

Storage upper-bandwidth.

type mu = int 

Upper-bandwidth.

type ml = int 

Lower-bandwidth.

type t = Sundials.RealArray2.t *
(smu * mu *
ml)

A band matrix accessible directly through a Bigarray.

The layout of these arrays are characterized by the storage upper bandwidth smu. Given an array a, the first dimension indexes the diagonals, with the main diagonal ($i = j$ ) at a.{smu, _}. The value in the ith row and jth column provided $\mathtt{i} \leq \mathtt{j} + \mathtt{ml}$ and $\mathtt{j} \leq \mathtt{i} + \mathtt{smu}$ is at a.{i - j + smu, j}.

Basic access

val make : smu * mu *
ml -> int -> float -> t

make (smu, mu, ml) n v returns an n by n band matrix with storage upper bandwidth smu, upper bandwidth sm, lower half-bandwidth ml, and all elements initialized to v.

If the result will not be LU factored then $\mathtt{smu} = \mathtt{mu}$ , otherwise $\mathtt{smu} = \min(\mathtt{n}-1, \mathtt{mu} + \mathtt{ml})$ . The extra space is used to store U after a call to Sundials_Matrix.ArrayBand.gbtrf.

val create : smu * mu *
ml -> int -> t

create smu ml n returns an uninitialized n by n band matrix with storage upper bandwidth smu and lower half-bandwidth ml.

val size : t -> int * int

m, n = size a returns the numbers of rows m and columns n of a.

NB: m and n are always equal for band matrices.

val dims : t ->
smu * mu *
ml

Returns the dimensions of an array band matrix.

val pp : Stdlib.Format.formatter -> t -> unit

Pretty-print a band matrix using the Format module.

val ppi : ?start:string ->
?stop:string ->
?sep:string ->
?indent:int ->
?itemsep:string ->
?empty:string ->
?item:(Stdlib.Format.formatter -> int -> int -> float -> unit) ->
unit -> Stdlib.Format.formatter -> t -> unit

Pretty-print an array band matrix using the Format module. The defaults are: start="[", stop="]", sep=";", indent=4, itemsep=" ", empty="           ~           " and item=fun f r c->Format.fprintf f "(%2d,%2d)=% -15e" r c (see fprintf). The indent argument specifies the indent for wrapped rows.

val get : t -> int -> int -> float

get a i j returns the value at row i and column j of a. Only rows and columns satisfying $\mathtt{i} \leq \mathtt{j} + \mathtt{ml}$ and $\mathtt{j} \leq \mathtt{i} + \mathtt{smu}$ are valid.

val set : t -> int -> int -> float -> unit

set a i j v sets the value at row i and column j of a to v. Only rows and columns satisfying $\mathtt{i} \leq \mathtt{j} + \mathtt{ml}$ and $\mathtt{j} \leq \mathtt{i} + \mathtt{smu}$ are valid.

val update : t -> int -> int -> (float -> float) -> unit

update a i j f sets the value at row i and column j of a to f v. Only rows and columns satisfying $\mathtt{i} \leq \mathtt{j} + \mathtt{ml}$ and $\mathtt{j} \leq \mathtt{i} + \mathtt{smu}$ are valid.

val unwrap : t -> Sundials.RealArray2.data

Direct access to the underlying storage array, which is accessed column first (unlike in Sundials_Matrix.ArrayBand.get).

Operations

val ops : (t, Sundials.RealArray.t)
Sundials_Matrix.matrix_ops

Operations on array-based band matrices.

val scale_add : float -> t -> t -> unit

scale_add c a b calculates $A = cA + B$.

NB: Unlike the Sundials_Matrix.Band.scale_add operation, this operation raises an exception if b has a greater bandwidth than a, i.e., it never resizes a.

val scale_addi : float -> t -> unit

scale_addi ml c A calculates $A = cA + I$.

val matvec : t ->
Sundials.RealArray.t -> Sundials.RealArray.t -> unit

The call matvec a x y computes the matrix-vector product $y = Ax$.

val set_to_zero : t -> unit

Fills the matrix with zeros.

val blit : src:t -> dst:t -> unit

blit ~src ~dst copies the contents of src into dst.

val space : t -> int * int

lrw, liw = space a returns the storage requirements of a as lrw realtype words and liw integer words.

Calculations

val add_identity : t -> unit

Increment a square matrix by the identity matrix.

val scale : float -> t -> unit

scale c a multiplies each element of the band matrix a by c.

val gbtrf : t -> Sundials.LintArray.t -> unit

gbtrf a p performs the LU factorization of a with partial pivoting according to p. The values in a are overwritten with those of the calculated L and U matrices. The diagonal belongs to U. The diagonal of L is all 1s. U may occupy elements up to bandwidth smu (rather than to mu).

val gbtrs : t ->
Sundials.LintArray.t -> Sundials.RealArray.t -> unit

gbtrs a p b finds the solution of ax = b using LU factorization. Both p and b must have the same number of rows as a.