Module Gsl_linalg


module Gsl_linalg: sig .. end
Simple linear algebra operations


Simple matrix multiplication


val matmult : a:Gsl_vectmat.mat ->
?transpa:bool ->
b:Gsl_vectmat.mat -> ?transpb:bool -> Gsl_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 : Gsl_vectmat.mat -> Gsl_permut.permut -> int
val _LU_solve : Gsl_vectmat.mat ->
Gsl_permut.permut -> b:Gsl_vectmat.vec -> x:Gsl_vectmat.vec -> unit
val _LU_svx : Gsl_vectmat.mat -> Gsl_permut.permut -> Gsl_vectmat.vec -> unit
val _LU_refine : a:Gsl_vectmat.mat ->
lu:Gsl_vectmat.mat ->
Gsl_permut.permut ->
b:Gsl_vectmat.vec -> x:Gsl_vectmat.vec -> res:Gsl_vectmat.vec -> unit
val _LU_invert : Gsl_vectmat.mat -> Gsl_permut.permut -> Gsl_vectmat.mat -> unit
val _LU_det : Gsl_vectmat.mat -> int -> float
val _LU_lndet : Gsl_vectmat.mat -> float
val _LU_sgndet : Gsl_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 ->
[< `A of float array * int * int
| `AA of float array array
| `M of Gsl_matrix.matrix
| `MF of Gsl_matrix_flat.matrix ] ->
Gsl_vectmat.mat * Gsl_permut.permut * int
val solve_LU : ?protect:bool ->
[< `A of float array * int * int
| `AA of float array array
| `M of Gsl_matrix.matrix
| `MF of Gsl_matrix_flat.matrix ] ->
[< `A of float array
| `V of Gsl_vector.vector
| `VF of Gsl_vector_flat.vector ] ->
float array
val det_LU : ?protect:bool ->
[< `A of float array * int * int
| `AA of float array array
| `M of Gsl_matrix.matrix
| `MF of Gsl_matrix_flat.matrix ] ->
float
val invert_LU : ?protect:bool ->
?result:Gsl_vectmat.mat ->
[< `A of float array * int * int
| `AA of float array array
| `M of Gsl_matrix.matrix
| `MF of Gsl_matrix_flat.matrix ] ->
Gsl_vectmat.mat

Complex LU decomposition


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

QR decomposition


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

QR Decomposition with Column Pivoting


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

Singular Value Decomposition


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

LQ decomposition


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

P^T L Q decomposition


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

Cholesky decomposition


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

Tridiagonal Decomposition of Real Symmetric Matrices


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

Tridiagonal Decomposition of Hermitian Matrices


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

Bidiagonalization


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

Householder solver


val _HH_solve : Gsl_vectmat.mat -> b:Gsl_vectmat.vec -> x:Gsl_vectmat.vec -> unit
val _HH_svx : Gsl_vectmat.mat -> Gsl_vectmat.vec -> unit
val solve_HH : ?protect:bool ->
[< `A of float array * int * int
| `AA of float array array
| `M of Gsl_matrix.matrix
| `MF of Gsl_matrix_flat.matrix ] ->
[< `A of float array
| `V of Gsl_vector.vector
| `VF of Gsl_vector_flat.vector ] ->
float array

Tridiagonal Systems


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

Exponential


val _exponential : Gsl_vectmat.mat -> Gsl_vectmat.mat -> Gsl_fun.mode -> unit
val exponential : ?mode:Gsl_fun.mode ->
[< `A of float array * int * int
| `M of Gsl_matrix.matrix
| `MF of Gsl_matrix_flat.matrix ] ->
[ `M of Gsl_matrix.matrix ]