Module Gsl.Linalg

Simple matrix multiplication

external matmult : a:Vectmat.mat -> ?⁠transpa:bool -> b:Vectmat.mat -> ?⁠transpb:bool -> Vectmat.mat -> unit = "ml_gsl_linalg_matmult_mod"

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

external _LU_decomp : Vectmat.mat -> Permut.permut -> int = "ml_gsl_linalg_LU_decomp"
external _LU_solve : Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_LU_solve"
external _LU_svx : Vectmat.mat -> Permut.permut -> Vectmat.vec -> unit = "ml_gsl_linalg_LU_svx"
external _LU_refine : a:Vectmat.mat -> lu:Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> res:Vectmat.vec -> unit = "ml_gsl_linalg_LU_refine_bc" "ml_gsl_linalg_LU_refine"
external _LU_invert : Vectmat.mat -> Permut.permut -> Vectmat.mat -> unit = "ml_gsl_linalg_LU_invert"
external _LU_det : Vectmat.mat -> int -> float = "ml_gsl_linalg_LU_det"
external _LU_lndet : Vectmat.mat -> float = "ml_gsl_linalg_LU_lndet"
external _LU_sgndet : Vectmat.mat -> int -> int = "ml_gsl_linalg_LU_sgndet"

Higher-level functions

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

external complex_LU_decomp : Vectmat.cmat -> Permut.permut -> int = "ml_gsl_linalg_complex_LU_decomp"
external complex_LU_solve : Vectmat.cmat -> Permut.permut -> b:Vectmat.cvec -> x:Vectmat.cvec -> unit = "ml_gsl_linalg_complex_LU_solve"
external complex_LU_svx : Vectmat.cmat -> Permut.permut -> Vectmat.cvec -> unit = "ml_gsl_linalg_complex_LU_svx"
external complex_LU_refine : a:Vectmat.cmat -> lu:Vectmat.cmat -> Permut.permut -> b:Vectmat.cvec -> x:Vectmat.cvec -> res:Vectmat.cvec -> unit = "ml_gsl_linalg_complex_LU_refine_bc" "ml_gsl_linalg_complex_LU_refine"
external complex_LU_invert : Vectmat.cmat -> Permut.permut -> Vectmat.cmat -> unit = "ml_gsl_linalg_complex_LU_invert"
external complex_LU_det : Vectmat.cmat -> int -> Gsl_complex.complex = "ml_gsl_linalg_complex_LU_det"
external complex_LU_lndet : Vectmat.cmat -> float = "ml_gsl_linalg_complex_LU_lndet"
external complex_LU_sgndet : Vectmat.cmat -> int -> Gsl_complex.complex = "ml_gsl_linalg_complex_LU_sgndet"

QR decomposition

external _QR_decomp : Vectmat.mat -> Vectmat.vec -> unit = "ml_gsl_linalg_QR_decomp"
external _QR_solve : Vectmat.mat -> Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QR_solve"
external _QR_svx : Vectmat.mat -> Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QR_svx"
external _QR_lssolve : Vectmat.mat -> Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> res:Vectmat.vec -> unit = "ml_gsl_linalg_QR_lssolve"
external _QR_QTvec : Vectmat.mat -> Vectmat.vec -> v:Vectmat.vec -> unit = "ml_gsl_linalg_QR_QTvec"
external _QR_Qvec : Vectmat.mat -> Vectmat.vec -> v:Vectmat.vec -> unit = "ml_gsl_linalg_QR_Qvec"
external _QR_Rsolve : Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QR_Rsolve"
external _QR_Rsvx : Vectmat.mat -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QR_Rsvx"
external _QR_unpack : Vectmat.mat -> tau:Vectmat.vec -> q:Vectmat.mat -> r:Vectmat.mat -> unit = "ml_gsl_linalg_QR_unpack"
external _QR_QRsolve : Vectmat.mat -> r:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QR_QRsolve"
external _QR_update : Vectmat.mat -> r:Vectmat.mat -> w:Vectmat.vec -> v:Vectmat.vec -> unit = "ml_gsl_linalg_QR_update"
external _R_solve : r:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_R_solve"

QR Decomposition with Column Pivoting

external _QRPT_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> norm:Vectmat.vec -> int = "ml_gsl_linalg_QRPT_decomp"
external _QRPT_decomp2 : a:Vectmat.mat -> q:Vectmat.mat -> r:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> norm:Vectmat.vec -> int = "ml_gsl_linalg_QRPT_decomp2_bc" "ml_gsl_linalg_QRPT_decomp2"
external _QRPT_solve : qr:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QRPT_solve"
external _QRPT_svx : qr:Vectmat.mat -> tau:Vectmat.vec -> p:Permut.permut -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QRPT_svx"
external _QRPT_QRsolve : q:Vectmat.mat -> r:Vectmat.mat -> p:Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QRPT_QRsolve"
external _QRPT_update : q:Vectmat.mat -> r:Vectmat.mat -> p:Permut.permut -> u:Vectmat.vec -> v:Vectmat.vec -> unit = "ml_gsl_linalg_QRPT_update"
external _QRPT_Rsolve : qr:Vectmat.mat -> p:Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QRPT_Rsolve"
external _QRPT_Rsvx : qr:Vectmat.mat -> p:Permut.permut -> x:Vectmat.vec -> unit = "ml_gsl_linalg_QRPT_Rsolve"

Singular Value Decomposition

external _SV_decomp : a:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> work:Vectmat.vec -> unit = "ml_gsl_linalg_SV_decomp"
external _SV_decomp_mod : a:Vectmat.mat -> x:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> work:Vectmat.vec -> unit = "ml_gsl_linalg_SV_decomp_mod"
external _SV_decomp_jacobi : a:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> unit = "ml_gsl_linalg_SV_decomp_jacobi"
external _SV_solve : u:Vectmat.mat -> v:Vectmat.mat -> s:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_SV_solve"

LQ decomposition

external _LQ_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_decomp"
external _LQ_solve_T : lq:Vectmat.mat -> tau:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_solve_T"
external _LQ_svx_T : lq:Vectmat.mat -> tau:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_svx_T"
external _LQ_lssolve_T : lq:Vectmat.mat -> tau:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> res:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_lssolve_T"
external _LQ_Lsolve_T : lq:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_Lsolve_T"
external _LQ_Lsvx_T : lq:Vectmat.mat -> x:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_Lsvx_T"
external _L_solve_T : l:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_L_solve_T"
external _LQ_vecQ : lq:Vectmat.mat -> tau:Vectmat.vec -> v:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_vecQ"
external _LQ_vecQT : lq:Vectmat.mat -> tau:Vectmat.vec -> v:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_vecQT"
external _LQ_unpack : lq:Vectmat.mat -> tau:Vectmat.vec -> q:Vectmat.mat -> l:Vectmat.mat -> unit = "ml_gsl_linalg_LQ_unpack"
external _LQ_update : q:Vectmat.mat -> r:Vectmat.mat -> v:Vectmat.vec -> w:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_update"
external _LQ_LQsolve : q:Vectmat.mat -> l:Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_LQ_LQsolve"

P^T L Q decomposition

external _PTLQ_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> norm:Vectmat.vec -> int = "ml_gsl_linalg_PTLQ_decomp"
external _PTLQ_decomp2 : a:Vectmat.mat -> q:Vectmat.mat -> r:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> norm:Vectmat.vec -> int = "ml_gsl_linalg_PTLQ_decomp2_bc" "ml_gsl_linalg_PTLQ_decomp2"
external _PTLQ_solve_T : qr:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_PTLQ_solve_T"
external _PTLQ_svx_T : lq:Vectmat.mat -> tau:Vectmat.vec -> Permut.permut -> x:Vectmat.vec -> unit = "ml_gsl_linalg_PTLQ_svx_T"
external _PTLQ_LQsolve_T : q:Vectmat.mat -> l:Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_PTLQ_LQsolve_T"
external _PTLQ_Lsolve_T : lq:Vectmat.mat -> Permut.permut -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_PTLQ_Lsolve_T"
external _PTLQ_Lsvx_T : lq:Vectmat.mat -> Permut.permut -> x:Vectmat.vec -> unit = "ml_gsl_linalg_PTLQ_Lsvx_T"
external _PTLQ_update : q:Vectmat.mat -> l:Vectmat.mat -> Permut.permut -> v:Vectmat.vec -> w:Vectmat.vec -> unit = "ml_gsl_linalg_PTLQ_update"

Cholesky decomposition

external cho_decomp : Vectmat.mat -> unit = "ml_gsl_linalg_cholesky_decomp"
external cho_solve : Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_cholesky_solve"
external cho_svx : Vectmat.mat -> Vectmat.vec -> unit = "ml_gsl_linalg_cholesky_svx"
external cho_decomp_unit : Vectmat.mat -> Vectmat.vec -> unit = "ml_gsl_linalg_cholesky_decomp_unit"

Tridiagonal Decomposition of Real Symmetric Matrices

external symmtd_decomp : a:Vectmat.mat -> tau:Vectmat.vec -> unit = "ml_gsl_linalg_symmtd_decomp"
external symmtd_unpack : a:Vectmat.mat -> tau:Vectmat.vec -> q:Vectmat.mat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit = "ml_gsl_linalg_symmtd_unpack"
external symmtd_unpack_T : a:Vectmat.mat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit = "ml_gsl_linalg_symmtd_unpack_T"

Tridiagonal Decomposition of Hermitian Matrices

external hermtd_decomp : a:Vectmat.cmat -> tau:Vectmat.cvec -> unit = "ml_gsl_linalg_hermtd_decomp"
external hermtd_unpack : a:Vectmat.cmat -> tau:Vectmat.cvec -> q:Vectmat.cmat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit = "ml_gsl_linalg_hermtd_unpack"
external hermtd_unpack_T : a:Vectmat.cmat -> diag:Vectmat.vec -> subdiag:Vectmat.vec -> unit = "ml_gsl_linalg_hermtd_unpack_T"

Bidiagonalization

external bidiag_decomp : a:Vectmat.mat -> tau_u:Vectmat.vec -> tau_v:Vectmat.vec -> unit = "ml_gsl_linalg_bidiag_decomp"
external 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 = "ml_gsl_linalg_bidiag_unpack_bc" "ml_gsl_linalg_bidiag_unpack"
external bidiag_unpack2 : a:Vectmat.mat -> tau_u:Vectmat.vec -> tau_v:Vectmat.vec -> v:Vectmat.mat -> unit = "ml_gsl_linalg_bidiag_unpack2"
external bidiag_unpack_B : a:Vectmat.mat -> diag:Vectmat.vec -> superdiag:Vectmat.vec -> unit = "ml_gsl_linalg_bidiag_unpack_B"

Householder solver

external _HH_solve : Vectmat.mat -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_HH_solve"
external _HH_svx : Vectmat.mat -> Vectmat.vec -> unit = "ml_gsl_linalg_HH_svx"
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

external solve_symm_tridiag : diag:Vectmat.vec -> offdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_solve_symm_tridiag"
external solve_tridiag : diag:Vectmat.vec -> abovediag:Vectmat.vec -> belowdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_solve_tridiag"
external solve_symm_cyc_tridiag : diag:Vectmat.vec -> offdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_solve_symm_cyc_tridiag"
external solve_cyc_tridiag : diag:Vectmat.vec -> abovediag:Vectmat.vec -> belowdiag:Vectmat.vec -> b:Vectmat.vec -> x:Vectmat.vec -> unit = "ml_gsl_linalg_solve_cyc_tridiag"

Exponential

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