Module Lacaml.C

Single precision complex BLAS and LAPACK functions.

This module Lacaml.C contains linear algebra routines for complex numbers (precision: complex32). It is recommended to use this module by writing

    open Lacaml.C

at the top of your file.

type prec = Bigarray.complex32_elt
type num_type = Complex.t
type vec = (Complex.t, Bigarray.complex32_elt, Bigarray.fortran_layout) Bigarray.Array1.t

Complex vectors (precision: complex32).

type rvec = (float, Bigarray.float32_elt, Bigarray.fortran_layout) Bigarray.Array1.t

Vectors of reals (precision: float32).

type mat = (Complex.t, Bigarray.complex32_elt, Bigarray.fortran_layout) Bigarray.Array2.t

Complex matrices (precision: complex32).

type trans3 = [
| `C
| `N
| `T
]

Transpose parameter (conjugate transposed, normal, or transposed).

val prec : (Complex.t, Bigarray.complex32_elt) Bigarray.kind

Precision for this submodule C. Allows to write precision independent code.

module Vec : sig ... end
module Mat : sig ... end
val pp_num : Format.formatter ‑> Complex.t ‑> unit

pp_num ppf el is equivalent to fprintf ppf "(%G, %Gi)" el.re el.im.

val pp_vec : (Complex.t, 'aIo.pp_vec

Pretty-printer for column vectors.

val pp_mat : (Complex.t, 'aIo.pp_mat

Pretty-printer for matrices.

BLAS-1 interface
val dotu : ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> vec ‑> Complex.t

dotu ?n ?ofsx ?incx x ?ofsy ?incy y see BLAS documentation!

val dotc : ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> vec ‑> Complex.t

dotc ?n ?ofsx ?incx x ?ofsy ?incy y see BLAS documentation!

LAPACK interface
val lansy_min_lwork : int ‑> Common.norm4 ‑> int

lansy_min_lwork m norm

val lansy : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠norm:Common.norm4 ‑> ?⁠work:rvec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> float

lansy ?n ?up ?norm ?work ?ar ?ac a see LAPACK documentation!

val gecon_min_lwork : int ‑> int

gecon_min_lwork n

val gecon_min_lrwork : int ‑> int

gecon_min_lrwork n

val gecon : ?⁠n:int ‑> ?⁠norm:Common.norm2 ‑> ?⁠anorm:float ‑> ?⁠work:vec ‑> ?⁠rwork:rvec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> float

gecon ?n ?norm ?anorm ?work ?rwork ?ar ?ac a

val sycon_min_lwork : int ‑> int

sycon_min_lwork n

val sycon : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠anorm:float ‑> ?⁠work:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> float

sycon ?n ?up ?ipiv ?anorm ?work ?ar ?ac a

val pocon_min_lwork : int ‑> int

pocon_min_lwork n

val pocon_min_lrwork : int ‑> int

pocon_min_lrwork n

val pocon : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠anorm:float ‑> ?⁠work:vec ‑> ?⁠rwork:rvec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> float

pocon ?n ?up ?anorm ?work ?rwork ?ar ?ac a

General Schur factorization
val gees : ?⁠n:int ‑> ?⁠jobvs:Common.schur_vectors ‑> ?⁠sort:Common.eigen_value_sort ‑> ?⁠w:vec ‑> ?⁠vsr:int ‑> ?⁠vsc:int ‑> ?⁠vs:mat ‑> ?⁠work:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> int * vec * mat

gees ?n ?jobvs ?sort ?w ?vsr ?vsc ?vs ?work ?ar ?ac a See gees-function for details about arguments.

General SVD routines
val gesvd_min_lwork : m:int ‑> n:int ‑> int

gesvd_min_lwork ~m ~n

val gesvd_lrwork : m:int ‑> n:int ‑> int

gesvd_lrwork m n

val gesvd_opt_lwork : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠jobu:Common.svd_job ‑> ?⁠jobvt:Common.svd_job ‑> ?⁠s:rvec ‑> ?⁠ur:int ‑> ?⁠uc:int ‑> ?⁠u:mat ‑> ?⁠vtr:int ‑> ?⁠vtc:int ‑> ?⁠vt:mat ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> int
val gesvd : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠jobu:Common.svd_job ‑> ?⁠jobvt:Common.svd_job ‑> ?⁠s:rvec ‑> ?⁠ur:int ‑> ?⁠uc:int ‑> ?⁠u:mat ‑> ?⁠vtr:int ‑> ?⁠vtc:int ‑> ?⁠vt:mat ‑> ?⁠work:vec ‑> ?⁠rwork:rvec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> rvec * mat * mat
General eigenvalue problem (simple drivers)
val geev_min_lwork : int ‑> int

geev_min_lwork n

val geev_min_lrwork : int ‑> int

geev_min_lrwork n

val geev_opt_lwork : ?⁠n:int ‑> ?⁠vlr:int ‑> ?⁠vlc:int ‑> ?⁠vl:mat option ‑> ?⁠vrr:int ‑> ?⁠vrc:int ‑> ?⁠vr:mat option ‑> ?⁠ofsw:int ‑> ?⁠w:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> int

geev ?work ?rwork ?n ?vlr ?vlc ?vl ?vrr ?vrc ?vr ?ofsw w ?ar ?ac a See geev-function for details about arguments.

val geev : ?⁠n:int ‑> ?⁠work:vec ‑> ?⁠rwork:vec ‑> ?⁠vlr:int ‑> ?⁠vlc:int ‑> ?⁠vl:mat option ‑> ?⁠vrr:int ‑> ?⁠vrc:int ‑> ?⁠vr:mat option ‑> ?⁠ofsw:int ‑> ?⁠w:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> mat * vec * mat

geev ?work ?rwork ?n ?vlr ?vlc ?vl ?vrr ?vrc ?vr ?ofsw w ?ar ?ac a

BLAS-1 interface
val swap : ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> vec ‑> unit

swap ?n ?ofsx ?incx x ?ofsy ?incy y see BLAS documentation!

val scal : ?⁠n:int ‑> Complex.t ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> unit

scal ?n alpha ?ofsx ?incx x see BLAS documentation!

val copy : ?⁠n:int ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> ?⁠y:vec ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> vec

copy ?n ?ofsy ?incy ?y ?ofsx ?incx x see BLAS documentation!

val nrm2 : ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> float

nrm2 ?n ?ofsx ?incx x see BLAS documentation!

val axpy : ?⁠alpha:Complex.t ‑> ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> vec ‑> unit

axpy ?alpha ?n ?ofsx ?incx x ?ofsy ?incy y see BLAS documentation!

val iamax : ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> int

iamax ?n ?ofsx ?incx x see BLAS documentation!

val amax : ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> Complex.t

amax ?n ?ofsx ?incx x

BLAS-2 interface
val gemv : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠beta:Complex.t ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> ?⁠y:vec ‑> ?⁠trans:trans3 ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> vec

gemv ?m ?n ?beta ?ofsy ?incy ?y ?trans ?alpha ?ar ?ac a ?ofsx ?incx x performs the operation y := alpha * op(a) * x + beta * y where op(a) = a or aᵀ according to the value of trans. See BLAS documentation for more information. BEWARE that the 1988 BLAS-2 specification mandates that this function has no effect when n=0 while the mathematically expected behavior is y ← beta * y.

val gbmv : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠beta:Complex.t ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> ?⁠y:vec ‑> ?⁠trans:trans3 ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> int ‑> int ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> vec

gbmv ?m ?n ?beta ?ofsy ?incy ?y ?trans ?alpha ?ar ?ac a kl ku ?ofsx ?incx x see BLAS documentation!

val symv : ?⁠n:int ‑> ?⁠beta:Complex.t ‑> ?⁠ofsy:int ‑> ?⁠incy:int ‑> ?⁠y:vec ‑> ?⁠up:bool ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> vec

symv ?n ?beta ?ofsy ?incy ?y ?up ?alpha ?ar ?ac a ?ofsx ?incx x see BLAS documentation!

val trmv : ?⁠n:int ‑> ?⁠trans:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> unit

trmv ?n ?trans ?diag ?up ?ar ?ac a ?ofsx ?incx x see BLAS documentation!

val trsv : ?⁠n:int ‑> ?⁠trans:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> unit

trsv ?n ?trans ?diag ?up ?ar ?ac a ?ofsx ?incx x see BLAS documentation!

val tpmv : ?⁠n:int ‑> ?⁠trans:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠up:bool ‑> ?⁠ofsap:int ‑> vec ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> unit

tpmv ?n ?trans ?diag ?up ?ofsap ap ?ofsx ?incx x see BLAS documentation!

val tpsv : ?⁠n:int ‑> ?⁠trans:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠up:bool ‑> ?⁠ofsap:int ‑> vec ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> unit

tpsv ?n ?trans ?diag ?up ?ofsap ap ?ofsx ?incx x see BLAS documentation!

BLAS-3 interface
val gemm : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠k:int ‑> ?⁠beta:Complex.t ‑> ?⁠cr:int ‑> ?⁠cc:int ‑> ?⁠c:mat ‑> ?⁠transa:trans3 ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠transb:trans3 ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> mat

gemm ?m ?n ?k ?beta ?cr ?cc ?c ?transa ?alpha ?ar ?ac a ?transb ?br ?bc b performs the operation c := alpha * op(a) * op(b) + beta * c where op(x) = x or xᵀ depending on transx. See BLAS documentation for more information.

val symm : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠side:Common.side ‑> ?⁠up:bool ‑> ?⁠beta:Complex.t ‑> ?⁠cr:int ‑> ?⁠cc:int ‑> ?⁠c:mat ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> mat

symm ?m ?n ?side ?up ?beta ?cr ?cc ?c ?alpha ?ar ?ac a ?br ?bc b see BLAS documentation!

val trmm : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠side:Common.side ‑> ?⁠up:bool ‑> ?⁠transa:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

trmm ?m ?n ?side ?up ?transa ?diag ?alpha ?ar ?ac a ?br ?bc b see BLAS documentation!

val trsm : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠side:Common.side ‑> ?⁠up:bool ‑> ?⁠transa:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

trsm ?m ?n ?side ?up ?transa ?diag ?alpha ?ar ?ac ~a ?br ?bc b see BLAS documentation!

val syrk : ?⁠n:int ‑> ?⁠k:int ‑> ?⁠up:bool ‑> ?⁠beta:Complex.t ‑> ?⁠cr:int ‑> ?⁠cc:int ‑> ?⁠c:mat ‑> ?⁠trans:Common.trans2 ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> mat

syrk ?n ?k ?up ?beta ?cr ?cc ?c ?trans ?alpha ?ar ?ac a see BLAS documentation!

val syr2k : ?⁠n:int ‑> ?⁠k:int ‑> ?⁠up:bool ‑> ?⁠beta:Complex.t ‑> ?⁠cr:int ‑> ?⁠cc:int ‑> ?⁠c:mat ‑> ?⁠trans:Common.trans2 ‑> ?⁠alpha:Complex.t ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> mat

syr2k ?n ?k ?up ?beta ?cr ?cc ?c ?trans ?alpha ?ar ?ac a ?br ?bc b see BLAS documentation!

LAPACK interface
Auxiliary routines
val lacpy : ?⁠uplo:[ `U | `L ] ‑> ?⁠patt:Common.Types.Mat.patt ‑> ?⁠m:int ‑> ?⁠n:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> ?⁠b:mat ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> mat

lacpy ?patt ?uplo ?m ?n ?br ?bc ?b ?ar ?ac a copy the (triangular) (sub-)matrix a (to an optional (sub-)matrix b) and return it. patt is more general than uplo and should be used in its place whenever strict BLAS conformance is not required. Only one of patt and uplo can be specified at a time.

val laswp : ?⁠n:int ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠k1:int ‑> ?⁠k2:int ‑> ?⁠incx:int ‑> Common.int32_vec ‑> unit

laswp ?n ?ar ?ac a ?k1 ?k2 ?incx ipiv swap rows of a according to ipiv. See LAPACK-documentation for details!

val lapmt : ?⁠forward:bool ‑> ?⁠m:int ‑> ?⁠n:int ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> Common.int32_vec ‑> unit

lapmt ?forward ?n ?m ?ar ?ac a k swap columns of a according to the permutations in k. See LAPACK-documentation for details!

val lassq : ?⁠n:int ‑> ?⁠scale:float ‑> ?⁠sumsq:float ‑> ?⁠ofsx:int ‑> ?⁠incx:int ‑> vec ‑> float * float

lassq ?n ?ofsx ?incx ?scale ?sumsq

val larnv : ?⁠idist:[ `Uniform0 | `Uniform1 | `Normal ] ‑> ?⁠iseed:Common.int32_vec ‑> ?⁠n:int ‑> ?⁠ofsx:int ‑> ?⁠x:vec ‑> unit ‑> vec

larnv ?idist ?iseed ?n ?ofsx ?x ()

val lange_min_lwork : int ‑> Common.norm4 ‑> int

lange_min_lwork m norm

val lange : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠norm:Common.norm4 ‑> ?⁠work:rvec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> float

lange ?m ?n ?norm ?work ?ar ?ac a

val lauum : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> unit

lauum ?n ?up ?ar ?ac a computes the product U * U**T or L**T * L, where the triangular factor U or L is stored in the upper or lower triangular part of the array a. The upper or lower part of a is overwritten.

Linear equations (computational routines)
val getrf : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> Common.int32_vec

getrf ?m ?n ?ipiv ?ar ?ac a computes an LU factorization of a general m-by-n matrix a using partial pivoting with row interchanges. See LAPACK documentation.

val getrs : ?⁠n:int ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠trans:trans3 ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

getrs ?n ?ipiv ?trans ?ar ?ac a ?nrhs ?br ?bc b solves a system of linear equations a * X = b or a' * X = b with a general n-by-n matrix a using the LU factorization computed by getrf. Note that matrix a will be passed to getrf if ipiv was not provided.

val getri_min_lwork : int ‑> int

getri_min_lwork n

val getri_opt_lwork : ?⁠n:int ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> int

getri_opt_lwork ?n ?ar ?ac a

val getri : ?⁠n:int ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠work:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> unit

getri ?n ?ipiv ?work ?ar ?ac a computes the inverse of a matrix using the LU factorization computed by getrf. Note that matrix a will be passed to getrf if ipiv was not provided.

val sytrf_min_lwork : unit ‑> int

sytrf_min_lwork ()

val sytrf_opt_lwork : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> int

sytrf_opt_lwork ?n ?up ?ar ?ac a

val sytrf : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠work:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> Common.int32_vec

sytrf ?n ?up ?ipiv ?work ?ar ?ac a computes the factorization of the real symmetric matrix a using the Bunch-Kaufman diagonal pivoting method.

val sytrs : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

sytrs ?n ?up ?ipiv ?ar ?ac a ?nrhs ?br ?bc b solves a system of linear equations a*X = b with a real symmetric matrix a using the factorization a = U*D*U**T or a = L*D*L**T computed by sytrf. Note that matrix a will be passed to sytrf if ipiv was not provided.

val sytri_min_lwork : int ‑> int

sytri_min_lwork n

val sytri : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠work:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> unit

sytri ?n ?up ?ipiv ?work ?ar ?ac a computes the inverse of the real symmetric indefinite matrix a using the factorization a = U*D*U**T or a = L*D*L**T computed by sytrf. Note that matrix a will be passed to sytrf if ipiv was not provided.

val potrf : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> unit

potrf ?n ?up ?ar ?ac a factorizes symmetric positive definite matrix a (or the designated submatrix) using Cholesky factorization.

val potrs : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

potrs ?n ?up ?ar ?ac a ?nrhs ?br ?bc b solves a system of linear equations a*X = b, where a is symmetric positive definite matrix, using the Cholesky factorization a = U**T*U or a = L*L**T computed by potrf.

val potri : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> unit

potri ?n ?up ?ar ?ac a computes the inverse of the real symmetric positive definite matrix a using the Cholesky factorization a = U**T*U or a = L*L**T computed by potrf.

val trtrs : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠trans:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

trtrs ?n ?up ?trans ?diag ?ar ?ac a ?nrhs ?br ?bc b solves a triangular system of the form a * X = b or a**T * X = n, where a is a triangular matrix of order n, and b is an n-by-nrhs matrix.

val tbtrs : ?⁠n:int ‑> ?⁠kd:int ‑> ?⁠up:bool ‑> ?⁠trans:trans3 ‑> ?⁠diag:Common.diag ‑> ?⁠abr:int ‑> ?⁠abc:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

tbtrs ?n ?kd ?up ?trans ?diag ?abr ?abc ab ?nrhs ?br ?bc b solves a triangular system of the form a * X = b or a**T * X = b, where a is a triangular band matrix of order n, and b is an n-by-nrhs matrix.

val trtri : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠diag:Common.diag ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> unit

trtri ?n ?up ?diag ?ar ?ac a computes the inverse of a real upper or lower triangular matrix a.

val geqrf_opt_lwork : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> int

geqrf_opt_lwork ?m ?n ?ar ?ac a

val geqrf_min_lwork : n:int ‑> int

geqrf_min_lwork ~n

val geqrf : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠work:vec ‑> ?⁠tau:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> vec

geqrf ?m ?n ?work ?tau ?ar ?ac a computes a QR factorization of a real m-by-n matrix a. See LAPACK documentation.

Linear equations (simple drivers)
val gesv : ?⁠n:int ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

gesv ?n ?ipiv ?ar ?ac a ?nrhs ?br ?bc b computes the solution to a real system of linear equations a * X = b, where a is an n-by-n matrix and X and b are n-by-nrhs matrices. The LU decomposition with partial pivoting and row interchanges is used to factor a as a = P * L * U, where P is a permutation matrix, L is unit lower triangular, and U is upper triangular. The factored form of a is then used to solve the system of equations a * X = b. On exit, b contains the solution matrix X.

val gbsv : ?⁠n:int ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠abr:int ‑> ?⁠abc:int ‑> mat ‑> int ‑> int ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

gbsv ?n ?ipiv ?abr ?abc ab kl ku ?nrhs ?br ?bc b computes the solution to a real system of linear equations a * X = b, where a is a band matrix of order n with kl subdiagonals and ku superdiagonals, and X and b are n-by-nrhs matrices. The LU decomposition with partial pivoting and row interchanges is used to factor a as a = L * U, where L is a product of permutation and unit lower triangular matrices with kl subdiagonals, and U is upper triangular with kl+ku superdiagonals. The factored form of a is then used to solve the system of equations a * X = b.

val gtsv : ?⁠n:int ‑> ?⁠ofsdl:int ‑> vec ‑> ?⁠ofsd:int ‑> vec ‑> ?⁠ofsdu:int ‑> vec ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

gtsv ?n ?ofsdl dl ?ofsd d ?ofsdu du ?nrhs ?br ?bc b solves the equation a * X = b where a is an n-by-n tridiagonal matrix, by Gaussian elimination with partial pivoting. Note that the equation A'*X = b may be solved by interchanging the order of the arguments du and dl.

val posv : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

posv ?n ?up ?ar ?ac a ?nrhs ?br ?bc b computes the solution to a real system of linear equations a * X = b, where a is an n-by-n symmetric positive definite matrix and X and b are n-by-nrhs matrices. The Cholesky decomposition is used to factor a as a = U**T * U, if up = true, or a = L * L**T, if up = false, where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of a is then used to solve the system of equations a * X = b.

val ppsv : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ofsap:int ‑> vec ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

ppsv ?n ?up ?ofsap ap ?nrhs ?br ?bc b computes the solution to the real system of linear equations a * X = b, where a is an n-by-n symmetric positive definite matrix stored in packed format and X and b are n-by-nrhs matrices. The Cholesky decomposition is used to factor a as a = U**T * U, if up = true, or a = L * L**T, if up = false, where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of a is then used to solve the system of equations a * X = b.

val pbsv : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠kd:int ‑> ?⁠abr:int ‑> ?⁠abc:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

pbsv ?n ?up ?kd ?abr ?abc ab ?nrhs ?br ?bc b computes the solution to a real system of linear equations a * X = b, where a is an n-by-n symmetric positive definite band matrix and X and b are n-by-nrhs matrices. The Cholesky decomposition is used to factor a as a = U**T * U, if up = true, or a = L * L**T, if up = false, where U is an upper triangular band matrix, and L is a lower triangular band matrix, with the same number of superdiagonals or subdiagonals as a. The factored form of a is then used to solve the system of equations a * X = b.

val ptsv : ?⁠n:int ‑> ?⁠ofsd:int ‑> vec ‑> ?⁠ofse:int ‑> vec ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

ptsv ?n ?ofsd d ?ofse e ?nrhs ?br ?bc b computes the solution to the real system of linear equations a*X = b, where a is an n-by-n symmetric positive definite tridiagonal matrix, and X and b are n-by-nrhs matrices. A is factored as a = L*D*L**T, and the factored form of a is then used to solve the system of equations.

val sysv_opt_lwork : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> int

sysv_opt_lwork ?n ?up ?ar ?ac a ?nrhs ?br ?bc b

val sysv : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠work:vec ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

sysv ?n ?up ?ipiv ?work ?ar ?ac a ?nrhs ?br ?bc b computes the solution to a real system of linear equations a * X = b, where a is an N-by-N symmetric matrix and X and b are n-by-nrhs matrices. The diagonal pivoting method is used to factor a as a = U * D * U**T, if up = true, or a = L * D * L**T, if up = false, where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The factored form of a is then used to solve the system of equations a * X = b.

val spsv : ?⁠n:int ‑> ?⁠up:bool ‑> ?⁠ipiv:Common.int32_vec ‑> ?⁠ofsap:int ‑> vec ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

spsv ?n ?up ?ipiv ?ofsap ap ?nrhs ?br ?bc b computes the solution to the real system of linear equations a * X = b, where a is an n-by-n symmetric matrix stored in packed format and X and b are n-by-nrhs matrices. The diagonal pivoting method is used to factor a as a = U * D * U**T, if up = true, or a = L * D * L**T, if up = false, where U (or L) is a product of permutation and unit upper (lower) triangular matrices, D is symmetric and block diagonal with 1-by-1 and 2-by-2 diagonal blocks. The factored form of a is then used to solve the system of equations a * X = b.

Least squares (simple drivers)
val gels_min_lwork : m:int ‑> n:int ‑> nrhs:int ‑> int

gels_min_lwork ~m ~n ~nrhs

val gels_opt_lwork : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠trans:Common.trans2 ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> int

gels_opt_lwork ?m ?n ?trans ?ar ?ac a ?nrhs ?br ?bc b

val gels : ?⁠m:int ‑> ?⁠n:int ‑> ?⁠work:vec ‑> ?⁠trans:Common.trans2 ‑> ?⁠ar:int ‑> ?⁠ac:int ‑> mat ‑> ?⁠nrhs:int ‑> ?⁠br:int ‑> ?⁠bc:int ‑> mat ‑> unit

gels ?m ?n ?work ?trans ?ar ?ac a ?nrhs ?br ?bc b see LAPACK documentation!