package gsl

  1. Overview
  2. Docs

Simple linear algebra operations

Simple matrix multiplication

val matmult : a:Vectmat.mat -> ?transpa:bool -> b:Vectmat.mat -> ?transpb:bool -> Vectmat.mat -> unit

matmult a ~transpa b ~transpb c stores in matrix c the product of matrices a and b. transpa or transpb allow transposition of either matrix, so it can compute a.b or Trans(a).b or a.Trans(b) or Trans(a).Trans(b) .

See also Gsl.Blas.gemm.

LU decomposition

Low-level functions
val _LU_decomp : Vectmat.mat -> Permut.permut -> int
val _LU_solve : Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _LU_svx : Vectmat.mat -> Permut.permut -> Vectmat.vec -> unit
val _LU_refine : a:Vectmat.mat -> lu:Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> res:Vectmat.vec -> unit
val _LU_invert : Vectmat.mat -> Permut.permut -> Vectmat.mat -> unit
val _LU_det : Vectmat.mat -> int -> float
val _LU_lndet : Vectmat.mat -> float
val _LU_sgndet : Vectmat.mat -> int -> int
Higher-level functions

With these, the arguments are protected (copied) and necessary intermediate datastructures are allocated;

val decomp_LU : ?protect:bool -> [< `M of Matrix.matrix | `MF of Matrix_flat.matrix | `A of float array * int * int | `AA of float array array ] -> Vectmat.mat * Permut.permut * int
val solve_LU : ?protect:bool -> [< `M of Matrix.matrix | `MF of Matrix_flat.matrix | `A of float array * int * int | `AA of float array array ] -> [< `A of float array | `VF of Vector_flat.vector | `V of Vector.vector ] -> float array
val det_LU : ?protect:bool -> [< `M of Matrix.matrix | `MF of Matrix_flat.matrix | `A of float array * int * int | `AA of float array array ] -> float
val invert_LU : ?protect:bool -> ?result:Vectmat.mat -> [< `M of Matrix.matrix | `MF of Matrix_flat.matrix | `A of float array * int * int | `AA of float array array ] -> Vectmat.mat

Complex LU decomposition

val complex_LU_decomp : Vectmat.cmat -> Permut.permut -> int
val complex_LU_solve : Vectmat.cmat -> Permut.permut -> b:Vectmat.cvec -> x:Vectmat.cvec -> unit
val complex_LU_svx : Vectmat.cmat -> Permut.permut -> Vectmat.cvec -> unit
val complex_LU_refine : a:Vectmat.cmat -> lu:Vectmat.cmat -> Permut.permut -> b:Vectmat.cvec -> x:Vectmat.cvec -> res:Vectmat.cvec -> unit
val complex_LU_invert : Vectmat.cmat -> Permut.permut -> Vectmat.cmat -> unit
val complex_LU_det : Vectmat.cmat -> int -> Gsl_complex.complex
val complex_LU_lndet : Vectmat.cmat -> float
val complex_LU_sgndet : Vectmat.cmat -> int -> Gsl_complex.complex

QR decomposition

val _QR_decomp : Vectmat.mat -> Vectmat.vec -> unit
val _QR_solve : Vectmat.mat -> Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _QR_svx : Vectmat.mat -> Vectmat.vec -> x:Vectmat.vec -> unit
val _QR_lssolve : Vectmat.mat -> Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> res:Vectmat.vec -> unit
val _QR_QTvec : Vectmat.mat -> Vectmat.vec -> v:Vectmat.vec -> unit
val _QR_Qvec : Vectmat.mat -> Vectmat.vec -> v:Vectmat.vec -> unit
val _QR_Rsolve : Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _QR_Rsvx : Vectmat.mat -> x:Vectmat.vec -> unit
val _QR_unpack : Vectmat.mat -> tau:Vectmat.vec -> q:Vectmat.mat -> r:Vectmat.mat -> unit
val _QR_QRsolve : Vectmat.mat -> r:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _QR_update : Vectmat.mat -> r:Vectmat.mat -> w:Vectmat.vec -> v:Vectmat.vec -> unit
val _R_solve : r:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit

QR Decomposition with Column Pivoting

val _QRPT_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> norm:Vectmat.vec -> int
val _QRPT_decomp2 : a:Vectmat.mat -> q:Vectmat.mat -> r:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> norm:Vectmat.vec -> int
val _QRPT_solve : qr:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _QRPT_svx : qr:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> x:Vectmat.vec -> unit
val _QRPT_QRsolve : q:Vectmat.mat -> r:Vectmat.mat -> p:Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _QRPT_update : q:Vectmat.mat -> r:Vectmat.mat -> p:Permut.permut -> u:Vectmat.vec -> v:Vectmat.vec -> unit
val _QRPT_Rsolve : qr:Vectmat.mat -> p:Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _QRPT_Rsvx : qr:Vectmat.mat -> p:Permut.permut -> x:Vectmat.vec -> unit

Singular Value Decomposition

val _SV_decomp : a:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> work:Vectmat.vec -> unit
val _SV_decomp_mod : a:Vectmat.mat -> x:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> work:Vectmat.vec -> unit
val _SV_decomp_jacobi : a:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> unit
val _SV_solve : u:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit

LQ decomposition

val _LQ_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> unit
val _LQ_solve_T : lq:Vectmat.mat -> tau:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _LQ_svx_T : lq:Vectmat.mat -> tau:Vectmat.vec -> x:Vectmat.vec -> unit
val _LQ_lssolve_T : lq:Vectmat.mat -> tau:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> res:Vectmat.vec -> unit
val _LQ_Lsolve_T : lq:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _LQ_Lsvx_T : lq:Vectmat.mat -> x:Vectmat.vec -> unit
val _L_solve_T : l:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _LQ_vecQ : lq:Vectmat.mat -> tau:Vectmat.vec -> v:Vectmat.vec -> unit
val _LQ_vecQT : lq:Vectmat.mat -> tau:Vectmat.vec -> v:Vectmat.vec -> unit
val _LQ_unpack : lq:Vectmat.mat -> tau:Vectmat.vec -> q:Vectmat.mat -> l:Vectmat.mat -> unit
val _LQ_update : q:Vectmat.mat -> r:Vectmat.mat -> v:Vectmat.vec -> w:Vectmat.vec -> unit
val _LQ_LQsolve : q:Vectmat.mat -> l:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit

P^T L Q decomposition

val _PTLQ_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> norm:Vectmat.vec -> int
val _PTLQ_decomp2 : a:Vectmat.mat -> q:Vectmat.mat -> r:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> norm:Vectmat.vec -> int
val _PTLQ_solve_T : qr:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _PTLQ_svx_T : lq:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> x:Vectmat.vec -> unit
val _PTLQ_LQsolve_T : q:Vectmat.mat -> l:Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _PTLQ_Lsolve_T : lq:Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _PTLQ_Lsvx_T : lq:Vectmat.mat -> Permut.permut -> x:Vectmat.vec -> unit
val _PTLQ_update : q:Vectmat.mat -> l:Vectmat.mat -> Permut.permut -> v:Vectmat.vec -> w:Vectmat.vec -> unit

Cholesky decomposition

val cho_decomp : Vectmat.mat -> unit
val cho_solve : Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val cho_svx : Vectmat.mat -> Vectmat.vec -> unit
val cho_decomp_unit : Vectmat.mat -> Vectmat.vec -> unit

Tridiagonal Decomposition of Real Symmetric Matrices

val symmtd_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> unit
val symmtd_unpack : a:Vectmat.mat -> tau:Vectmat.vec -> q:Vectmat.mat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit
val symmtd_unpack_T : a:Vectmat.mat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit

Tridiagonal Decomposition of Hermitian Matrices

val hermtd_decomp : a:Vectmat.cmat -> tau:Vectmat.cvec -> unit
val hermtd_unpack : a:Vectmat.cmat -> tau:Vectmat.cvec -> q:Vectmat.cmat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit
val hermtd_unpack_T : a:Vectmat.cmat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit

Bidiagonalization

val bidiag_decomp : a:Vectmat.mat -> tau_u:Vectmat.vec -> tau_v:Vectmat.vec -> unit
val bidiag_unpack : a:Vectmat.mat -> tau_u:Vectmat.vec -> u:Vectmat.mat -> tau_v:Vectmat.vec -> v:Vectmat.mat -> diag:Vectmat.vec -> superdiag:Vectmat.vec -> unit
val bidiag_unpack2 : a:Vectmat.mat -> tau_u:Vectmat.vec -> tau_v:Vectmat.vec -> v:Vectmat.mat -> unit
val bidiag_unpack_B : a:Vectmat.mat -> diag:Vectmat.vec -> superdiag:Vectmat.vec -> unit

Householder solver

val _HH_solve : Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val _HH_svx : Vectmat.mat -> Vectmat.vec -> unit
val solve_HH : ?protect:bool -> [< `M of Matrix.matrix | `MF of Matrix_flat.matrix | `A of float array * int * int | `AA of float array array ] -> [< `A of float array | `VF of Vector_flat.vector | `V of Vector.vector ] -> float array

Tridiagonal Systems

val solve_symm_tridiag : diag:Vectmat.vec -> offdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val solve_tridiag : diag:Vectmat.vec -> abovediag:Vectmat.vec -> belowdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val solve_symm_cyc_tridiag : diag:Vectmat.vec -> offdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit
val solve_cyc_tridiag : diag:Vectmat.vec -> abovediag:Vectmat.vec -> belowdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit

Exponential

val _exponential : Vectmat.mat -> Vectmat.mat -> Fun.mode -> unit
val exponential : ?mode:Fun.mode -> [< `M of Matrix.matrix | `MF of Matrix_flat.matrix | `A of float array * int * int ] -> [ `M of Matrix.matrix ]
OCaml

Innovation. Community. Security.