Math::Libgsl::LinearAlgebra

An interface to libgsl, the Gnu Scientific Library - Linear Algebra.

NAME

Math::Libgsl::LinearAlgebra - An interface to libgsl, the Gnu Scientific Library - Linear Algebra.

SYNOPSIS

use Math::Libgsl::LinearAlgebra;

DESCRIPTION

Math::Libgsl::LinearAlgebra is an interface to the linear algebra functions of libgsl, the GNU Scientific Library. This package provides both the low-level interface to the C library (Raw) and a more comfortable interface layer for the Raku programmer.

This module provides functions for Num and Complex data types.

Num

LU-decomp(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function factorizes the matrix A into the LU decomposition PA = LU. The factorization is done in place, so on output the diagonal and upper triangular (or trapezoidal) part of the input matrix A contain the matrix U. The lower triangular (or trapezoidal) part of the input matrix (excluding the diagonal) contains L. The diagonal elements of L are unity, and are not stored. The return value is a List: the sign of the permutation and a permutation object, which encodes the permutation matrix P. In case of error a failure object is returned.

LU-solve(Math::Libgsl::Matrix LU.matrix.size2, Math::Libgsl::Permutation LU.matrix.size1, Math::Libgsl::Vector LU.matrix.size1 --> Math::Libgsl::Vector)

This function solves the square system Ax = b using the LU decomposition of A into (LU, p) given by the output of LU-decomp. In case of error a failure object is returned.

LU-svx(Math::Libgsl::Matrix LU.matrix.size2, Math::Libgsl::Permutation LU.matrix.size1, Math::Libgsl::Vector LU.matrix.size1 --> Int)

This function solves the square system Ax = b in-place using the precomputed LU decomposition of A into (LU, p). On input $x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-refine(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Matrix LU.matrix.size1 == A.matrix.size1 == p where _.size == b where _.vector.size == x where _.vector.size == $LU.matrix.size1 --> Int)

This function applies an iterative improvement to x, the solution of Ax = b, from the precomputed LU decomposition of A into (LU, p). This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-invert(Math::Libgsl::Matrix p where *.size == $LU.matrix.size1 --> Math::Libgsl::Matrix)

This function computes the inverse of a matrix A from its LU decomposition (LU, p), returning the matrix inverse. In case of error a failure object is returned.

LU-det(Math::Libgsl::Matrix signum where * ~~ -1|1 --> Num)

This function computes the determinant of a matrix A from its LU decomposition, signum. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-lndet(Math::Libgsl::Matrix $LU --> Num)

This function computes the determinant the logarithm of the absolute value of the determinant of a matrix A, ln |det(A)|, from its LU decomposition, $LU. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-sgndet(Math::Libgsl::Matrix signum where * ~~ -1|1 --> Int)

This function computes the sign or phase factor of the determinant of a matrix A, det(A)/|det(A)| from its LU decomposition, $LU. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QR-decomp(Math::Libgsl::Matrix $A --> Math::Libgsl::Vector)

This function factorizes the M-by-N matrix $A into the QR decomposition A = QR. On output the diagonal and upper triangular part of the input matrix contain the matrix R. The returned vector and the columns of the lower triangular part of the matrix A contain the Householder coefficients and Householder vectors which encode the orthogonal matrix Q. In case of error a failure object is returned.

QR-solve(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1 --> Math::Libgsl::Vector)

This function solves the square system Ax = b using the QR decomposition of A held in (tau) which must have been computed previously with QR-decomp(). In case of error a failure object is returned.

QR-svx(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1 --> Int)

This function solves the square system Ax = b in-place using the QR decomposition of A held in (tau) which must have been computed previously by QR-decomp(). On input $x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QR-lssolve(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1 --> List)

This function finds the least squares solution to the overdetermined system Ax = b where the matrix A has more rows than columns. The least squares solution minimizes the Euclidean norm of the residual, ||Ax − b||. The routine requires as input the QR decomposition of A into (tau) given by QR-decomp(). The function returns a List of two Math::Libgsl::Vector objects: the solution x and the residual. In case of error a failure object is returned.

QR-QTvec(Math::Libgsl::Matrix tau where *.vector.size == min(QR.matrix.size2), Math::Libgsl::Vector QR.matrix.size1 --> Int)

These function applies the matrix T(Q) encoded in the decomposition (QR, tau) to the vector v. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix T(Q). This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QR-Qvec(Math::Libgsl::Matrix tau where *.vector.size == min(QR.matrix.size2), Math::Libgsl::Vector QR.matrix.size1 --> Int)

This function applies the matrix Q encoded in the decomposition (tau) to the vector v. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix Q. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QR-QTmat(Math::Libgsl::Matrix tau where *.vector.size == min(QR.matrix.size2), Math::Libgsl::Matrix QR.matrix.size1 --> Int)

This function applies the matrix T(Q) encoded in the decomposition (tau) to the M-by-K matrix B. The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix T(Q). This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QR-Rsolve(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1 --> Math::Libgsl::Vector)

This function solves the triangular system Rx = b and returns the Math::Libgsl::Vector object $x. In case of error a failure object is returned.

QR-Rsvx(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size2 --> Int)

This function solves the triangular system Rx = b for x in-place. On input $x should contain the right-hand side b and is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QR-unpack(Math::Libgsl::Matrix tau where *.vector.size == min(QR.matrix.size2) --> List)

This function unpacks the encoded QR decomposition (tau) into the matrices Q and R. The function returns a List of two Math::Libgsl::Matrix objects: R which is M-by-N. In case of error a failure object is returned.

QR-QRsolve(Math::Libgsl::Matrix R where { R.matrix.size2 && R.matrix.size1 }, Math::Libgsl::Vector R.matrix.size1 --> Math::Libgsl::Vector)

This function solves the system Rx = T(Q) b for x. It can be used when the QR decomposition of a matrix is available in unpacked form as (R). The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

QR-update(Math::Libgsl::Matrix R where { R.matrix.size1 && R.matrix.size1 }, Math::Libgsl::Vector R.matrix.size1, Math::Libgsl::Vector R.matrix.size2 --> Int)

This function performs a rank-1 update wT(v) of the QR decomposition (R). The update is given by Q'R' = Q(R+wT(v)) where the output matrices R are also orthogonal and right triangular. Note that $w is destroyed by the update. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

R-solve(Math::Libgsl::Matrix R.matrix.size2, Math::Libgsl::Vector R.matrix.size1 --> Math::Libgsl::Vector)

This function solves the triangular system Rx = b for the N-by-N matrix x. In case of error a failure object is returned.

R-svx(Math::Libgsl::Matrix R.matrix.size2, Math::Libgsl::Vector R.matrix.size2 --> Int)

This function solves the triangular system Rx = b in-place. On input $x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QRPT-decomp(Math::Libgsl::Matrix $A --> List)

This function factorizes the M-by-N matrix tau, the Math::Libgsl::Permutation signum. In case of error a failure object is returned.

QRPT-decomp2(Math::Libgsl::Matrix $A --> List)

This function factorizes the matrix A itself. The function returns a List: the Math::Libgsl::Matrix R, the Math::Libgsl::Permutation signum. In case of error a failure object is returned.

QRPT-solve(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1, p where _.p.size == b where _.vector.size == $QR.matrix.size1 --> Math::Libgsl::Vector)

This function solves the square system Ax = b using the QRT(P) decomposition of A held in (tau, x. In case of error a failure object is returned.

QRPT-svx(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1, p where _.p.size == x where _.vector.size == $QR.matrix.size2 --> Int)

This function solves the square system Ax = b in-place using the QRT(P) decomposition of A held in (tau, x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QRPT-lssolve(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1, p where _.p.size == b where _.vector.size == $QR.matrix.size1 --> List)

This function finds the least squares solution to the overdetermined system Ax = b where the matrix A has more rows than columns and is assumed to have full rank. The least squares solution minimizes the Euclidean norm of the residual, ||b − Ax||. The routine requires as input the QR decomposition of A into (tau, $p) given by QRPT-decomp. The function returns a List of two Math::Libgsl::Vector objects: the solution x and the residual. In case of error a failure object is returned.

QRPT-lssolve2(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1, p where _.p.size == b where _.vector.size == rank where 0 < * ≤ $QR.matrix.size2 --> List)

This function finds the least squares solution to the overdetermined system Ax = b where the matrix A has more rows than columns and has rank given by the input rank. If the user does not know the rank of A, it may be estimated by calling QRPT-rank. The routine requires as input the QR decomposition of A into (tau, $p) given by QRPT-decomp. The function returns a List of two Math::Libgsl::Vector objects: the solution x and the residual. In case of error a failure object is returned.

QRPT-QRsolve(Math::Libgsl::Matrix Q.matrix.size2, Math::Libgsl::Matrix R.matrix.size1 == R.matrix.size1 == p where _.p.size == b where _.vector.size == $Q.matrix.size1 --> Math::Libgsl::Vector)

This function solves the square system RT(P) x = T(Q) b for x. It can be used when the QR decomposition of a matrix is available in unpacked form as (R). The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

QRPT-update(Math::Libgsl::Matrix R where { R.matrix.size1 && R.matrix.size1 }, Math::Libgsl::Permutation R.matrix.size1, Math::Libgsl::Vector R.matrix.size1, Math::Libgsl::Vector R.matrix.size2 --> Int)

This function performs a rank-1 update wT(v) of the QRT(P) decomposition (R, w is destroyed by the update. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QRPT-Rsolve(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Permutation QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size1 --> Math::Libgsl::Vector)

This function solves the triangular system RT(P) x = b for the N-by-N matrix R contained in x. In case of error a failure object is returned.

QRPT-Rsvx(Math::Libgsl::Matrix QR.matrix.size2, Math::Libgsl::Permutation QR.matrix.size2, Math::Libgsl::Vector QR.matrix.size2 --> Int)

This function solves the triangular system RT(P) x = b in-place for the N-by-N matrix R contained in x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

QRPT-rank(Math::Libgsl::Matrix tolerance --> Int)

This function returns the rank of the triangular matrix R contained in $QR.

QRPT-rcond(Math::Libgsl::Matrix QR.matrix.size2 --> Num)

This function returns the reciprocal condition number (using the 1-norm) of the R factor, stored in the upper triangle of $QR. In case of error a failure object is returned.

LQ-decomp(Math::Libgsl::Matrix $A --> Math::Libgsl::Vector)

This function factorizes the M-by-N matrix tau. The vector A contain the Householder coefficients and Householder vectors which encode the orthogonal matrix Q. In case of error a failure object is returned.

LQ-solve-T(Math::Libgsl::Matrix LQ.matrix.size2, Math::Libgsl::Vector LQ.matrix.size1, b where *.vector.size == $LQ.matrix.size2 --> Math::Libgsl::Vector)

This function finds the solution to the system Ax = b. The routine requires as input the LQ decomposition of A into (tau) given by LQ-decomp. The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

LQ-svx-T(Math::Libgsl::Matrix LQ.matrix.size2, Math::Libgsl::Vector LQ.matrix.size1, $LQ.matrix.size2), --> Math::Libgsl::Vector)

This function finds the solution to the system Ax = b. The routine requires as input the LQ decomposition of A into (tau) given by LQ-decomp. On input $x should contain the right-hand side b, which is replaced by the solution on output. In case of error a failure object is returned.

LQ-lssolve-T(Math::Libgsl::Matrix LQ.matrix.size2, Math::Libgsl::Vector LQ.matrix.size1, b where *.vector.size == $LQ.matrix.size1 --> List)

This function finds the minimum norm least squares solution to the underdetermined system Ax = b, where the M-by-N matrix A has M ≤ N. The routine requires as input the LQ decomposition of A into (tau) given by LQ-decomp. The function returns a List of two Math::Libgsl::Vector objects: the solution x and the residual. In case of error a failure object is returned.

LQ-Lsolve-T(Math::Libgsl::Matrix LQ.matrix.size2, Math::Libgsl::Vector LQ.matrix.size1 --> Math::Libgsl::Vector)

The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

LQ-Lsvx-T(Math::Libgsl::Matrix LQ.matrix.size2 --> Math::Libgsl::Vector)

On input $x should contain the right-hand side b, which is replaced by the solution on output. In case of error a failure object is returned.

L-solve-T(Math::Libgsl::Matrix L.matrix.size2, Math::Libgsl::Vector L.matrix.size2 --> Math::Libgsl::Vector)

The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

LQ-vecQ(Math::Libgsl::Matrix tau where *.vector.size == min(LQ.matrix.size2), Math::Libgsl::Vector LQ.matrix.size1 --> Int)

This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LQ-vecQT(Math::Libgsl::Matrix tau where *.vector.size == min(LQ.matrix.size2), Math::Libgsl::Vector LQ.matrix.size1 --> Int)

This function applies T(Q) to the vector v, storing the result T(Q) v in $v on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LQ-unpack(Math::Libgsl::Matrix tau where *.vector.size == min(LQ.matrix.size2) --> List)

This function unpacks the encoded LQ decomposition (tau). The function outputs a List: the Math::Libgsl::Matrix L. In case of error a failure object is returned.

LQ-update(Math::Libgsl::Matrix L where { Q.matrix.size1 && Q.matrix.size2 }, Math::Libgsl::Vector L.matrix.size1, Math::Libgsl::Vector L.matrix.size2 --> Int)

This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LQ-LQsolve(Math::Libgsl::Matrix L where { L.matrix.size2 && L.matrix.size2 }, Math::Libgsl::Vector L.matrix.size2 --> Math::Libgsl::Vector)

In case of error a failure object is returned.

COD-decomp(Math::Libgsl::Matrix $A --> List)

COD-decomp-e(Math::Libgsl::Matrix tolerance --> List)

These functions factor the M-by-N matrix A. The matrices Q and Z are encoded in packed storage in A. In case of error a failure object is returned.

COD-lssolve(Math::Libgsl::Matrix QRZT.matrix.size2, Math::Libgsl::Vector QRZT.matrix.size1, tau-Z where *.vector.size == min(QRZT.matrix.size2), Math::Libgsl::Permutation QRZT.matrix.size2, Int QRZT.matrix.size1, b where *.vector.size == $QRZT.matrix.size1 --> List)

This function finds the unique minimum norm least squares solution to the overdetermined system Ax = b where the matrix A into (tau_Q, p, $rank) given by COD-decomp. The function outputs a List: a Math::Libgsl::Vector object which is the solution x, and a Math::Libgsl::Vector object which stores the residual b − Ax. In case of error a failure object is returned.

COD-lssolve2(Math::Libgsl::Matrix QRZT.matrix.size2, Math::Libgsl::Vector QRZT.matrix.size1, tau-Z where *.vector.size == min(QRZT.matrix.size2), Math::Libgsl::Permutation QRZT.matrix.size2, Int QRZT.matrix.size1, b where *.vector.size == lambda --> List)

This function finds the solution to the regularized least squares problem in Tikhonov standard form, minₓ||b − Ax||² + λ²||x||². The routine requires as input the QRZT decomposition of A into (tau_Q, p, lambda. The function outputs a List: a Math::Libgsl::Vector object which is the solution x, and a Math::Libgsl::Vector object which stores the residual b − Ax. In case of error a failure object is returned.

COD-unpack(Math::Libgsl::Matrix tau-Q where *.vector.size == min(QRZT.matrix.size2), Math::Libgsl::Vector QRZT.matrix.size1, rank where * ≤ min(QRZT.matrix.size2) --> List)

This function unpacks the encoded QRZT decomposition (tau_Q, rank). The function returns a List of three Math::Libgsl::Matrix objects: R, $Z. In case of error a failure object is returned.

COD-matZ(Math::Libgsl::Matrix tau-Z where *.vector.size == min(QRZT.matrix.size2), Math::Libgsl::Matrix QRZT.matrix.size2, Int $rank --> Int)

This function multiplies the input matrix QRZT, rank). This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

SV-decomp(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function factorizes the M-by-N matrix A is replaced by U. The function returns a List: a Math::Libgsl::Matrix object which contains the elements of V in untransposed form, and a Math::Libgsl::Vector object which contains the diagonal elements of the singular value matrix S. The singular values are non-negative and form a non-increasing sequence from S₁ to Sₙ. In case of error a failure object is returned.

SV-decomp-mod(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function computes the SVD using the modified Golub-Reinsch algorithm, which is faster for M ≫ N. The function returns a List: a Math::Libgsl::Matrix object which contains the elements of V in untransposed form, and a Math::Libgsl::Vector object which contains the diagonal elements of the singular value matrix S. The singular values are non-negative and form a non-increasing sequence from S₁ to Sₙ. In case of error a failure object is returned.

SV-decomp-jacobi(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function computes the SVD of the M-by-N matrix A using one-sided Jacobi orthogonalization for M ≥ N. The function returns a List: a Math::Libgsl::Matrix object which contains the elements of V in untransposed form, and a Math::Libgsl::Vector object which contains the diagonal elements of the singular value matrix S. The singular values are non-negative and form a non-increasing sequence from S₁ to Sₙ. In case of error a failure object is returned.

SV-solve(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Matrix V.matrix.size1 == V.matrix.size1 == S where _.vector.size == b where *.vector.size == $A.matrix.size1 --> Math::Libgsl::Vector)

This function solves the system Ax = b using the singular value decomposition (U, S, V) of A which must have been computed previously with COD-decomp. Only non-zero singular values are used in computing the solution. The parts of the solution corresponding to singular values of zero are ignored. Other singular values can be edited out by setting them to zero before calling this function. The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

SV-leverage(Math::Libgsl::Matrix $A --> Math::Libgsl::Vector)

This function computes the statistical leverage values hᵢ of a matrix A using its singular value decomposition (U, S, V) previously computed with COD-decomp. The function returns a Math::Libgsl::Vector object which stores the diagonal values of the matrix A(T(A)A)⁻¹T(A) and depend only on the matrix U which is the input to this function. In case of error a failure object is returned.

cholesky-decomp1(Math::Libgsl::Matrix A.matrix.size2 --> Int)

This function factorizes the symmetric, positive-definite square matrix A are used (the upper triangular part is ignored). On output the diagonal and lower triangular part of the input matrix $A contain the matrix L, while the upper triangular part contains the original matrix. If the matrix is not positive-definite then the decomposition will fail, returning the error code GSL_EDOM. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-solve(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size1 --> Math::Libgsl::Vector)

These functions solve the system Ax = b using the Cholesky decomposition of x. In case of error a failure object is returned.

cholesky-svx(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size1 --> Int)

This function solves the system Ax = b in-place using the Cholesky decomposition of x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-solve-mat(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Matrix A.matrix.size1 --> Math::Libgsl::Matrix)

In case of error a failure object is returned.

cholesky-svx-mat(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Matrix A.matrix.size1 --> Int)

This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-invert(Math::Libgsl::Matrix A.matrix.size2 --> Int)

These functions compute the inverse of a matrix from its Cholesky decomposition A. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-decomp2(Math::Libgsl::Matrix A.matrix.size2 --> Math::Libgsl::Vector)

This function calculates a diagonal scaling transformation S for the symmetric, positive-definite square matrix A are used (the upper triangular part is ignored). On output the diagonal and lower triangular part of the input matrix $A contain the matrix L, while the upper triangular part of the input matrix is overwritten with T(L) (the diagonal terms being identical for both L and L T). If the matrix is not positive-definite then the decomposition will fail, returning the error code GSL_EDOM. The function returns a Math::Libgsl::Vector object which stores the diagonal scale factors. In case of error a failure object is returned.

cholesky-solve2(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size1, Math::Libgsl::Vector A.matrix.size1 --> Math::Libgsl::Vector)

This function solves the system (SAS)(S⁻¹ x) = Sb using the Cholesky decomposition of SAS held in the matrix x. In case of error a failure object is returned.

cholesky-svx2(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size2, Math::Libgsl::Vector A.matrix.size2 --> Int)

This function solves the system (SAS)(S⁻¹ x) = Sb using the Cholesky decomposition of SAS held in the matrix x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-decomp-unit(Math::Libgsl::Matrix A.matrix.size2 --> Math::Libgsl::Vector)

The function returns the Math::Libgsl::Vector $d. In case of error a failure object is returned.

cholesky-scale(Math::Libgsl::Matrix A.matrix.size2 --> Math::Libgsl::Vector)

This function calculates a diagonal scaling transformation of the symmetric, positive definite matrix $A, such that SAS has a condition number within a factor of N of the matrix of smallest possible condition number over all possible diagonal scalings. The function outputs a Math::Libgsl::Vector object which contains the scale factors. In case of error a failure object is returned.

cholesky-scale-apply(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size2 --> Int)

This function applies the scaling transformation S to the matrix A is replaced by SAS. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-rcond(Math::Libgsl::Matrix A.matrix.size2 --> Num)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric positive definite matrix A, using its Cholesky decomposition provided in $A. The function returns a Math::Libgsl::Vector object which stores the reciprocal condition number estimate, defined as 1/(||A||₁ · ||A⁻¹||₁). In case of error a failure object is returned.

pcholesky-decomp(Math::Libgsl::Matrix A.matrix.size2 --> Math::Libgsl::Permutation)

This function factors the symmetric, positive-definite square matrix A are used to construct the factorization. On output the diagonal of the input matrix A contains the matrix L. Since L has ones on its diagonal these do not need to be explicitely stored. The upper triangular portion of $A is unmodified. The function returns the Math::Libgsl::Permutation object which stores the permutation matrix P. In case of error a failure object is returned.

pcholesky-solve(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1 --> Math::Libgsl::Vector)

This function solves the system Ax = b using the Pivoted Cholesky decomposition of A held in the matrix p which must have been previously computed by pcholesky-decomp. The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

pcholesky-svx(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1 --> Int)

This function solves the system Ax = b using the Pivoted Cholesky decomposition of A held in the matrix p which must have been previously computed by pcholesky-decomp. On input, $x contains the right hand side vector b which is replaced by the solution vector on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

pcholesky-decomp2(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function computes the pivoted Cholesky factorization of the matrix SAS, where the input matrix A as much as possible. On input, the values from the diagonal and lower-triangular part of the matrix A stores the diagonal elements of D, and the lower triangular portion of A is unmodified. The function returns a List: the permutation matrix P is stored in a Math::Libgsl::Permutation object, the diagonal scaling transformation is stored in a Math::Libgsl::Vector object. In case of error a failure object is returned.

pcholesky-solve2(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1 --> Math::Libgsl::Vector)

This function solves the system (SAS)(S⁻¹x) = Sb using the Pivoted Cholesky decomposition of SAS held in the matrix p, and vector x. In case of error a failure object is returned.

pcholesky-svx2(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1 --> Int)

This function solves the system (SAS)(S⁻¹x) = Sb using the Pivoted Cholesky decomposition of SAS held in the matrix p, and vector x contains the right hand side vector b which is replaced by the solution vector on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

pcholesky-invert(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1 --> Math::Libgsl::Matrix)

This function computes the inverse of the matrix A, using the Pivoted Cholesky decomposition stored in p. The function returns the Math::Libgsl::Matrix A⁻¹. In case of error a failure object is returned.

pcholesky-rcond(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1 --> Num)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric positive definite matrix A, using its pivoted Cholesky decomposition provided in $LDLT. The function returns the reciprocal condition number estimate, defined as 1/(||A||₁ · ||A⁻¹||₁). This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

mcholesky-decomp(Math::Libgsl::Matrix A.matrix.size2, Bool :$perturbation = True --> List)

This function factors the symmetric, indefinite square matrix A are used to construct the factorization. On output the diagonal of the input matrix A contains the matrix L. Since L has ones on its diagonal these do not need to be explicitely stored. The upper triangular portion of $A is unmodified. The function returns a List: the permutation matrix P, stored in a Math::Libgsl::Permutation object and the diagonal perturbation matrix, stored in a Math::Libgsl::Vector object. In case of error a failure object is returned.

mcholesky-solve(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1 --> Math::Libgsl::Vector)

This function solves the perturbed system (A + E)x = b using the Cholesky decomposition of A + E held in the matrix p which must have been previously computed by mcholesky-decomp. The function returns the Math::Libgsl::Vector $x. In case of error a failure object is returned.

mcholesky-svx(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1, Math::Libgsl::Vector LDLT.matrix.size1 --> Int)

This function solves the perturbed system (A + E)x = b using the Cholesky decomposition of A + E held in the matrix p which must have been previously computed by mcholesky-decomp. On input, $x contains the right hand side vector b which is replaced by the solution vector on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

mcholesky-rcond(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1 --> Num)

This function estimates the reciprocal condition number (using the 1-norm) of the perturbed matrix A+E, using its pivoted Cholesky decomposition provided in $LDLT. The function returns the reciprocal condition number estimate, defined as 1/(||A + E||₁ · ||(A + E)⁻¹||₁). In case of error a failure object is returned.

mcholesky-invert(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Permutation LDLT.matrix.size1 --> Math::Libgsl::Matrix)

The function returns the inverted matrix. In case of error a failure object is returned.

ldlt-decomp(Math::Libgsl::Matrix A.matrix.size2 --> Int)

This function factorizes the symmetric, non-singular square matrix A are used. The upper triangle of A contains the matrix D and the lower triangle of $A contains the unit lower triangular matrix L. The matrix 1-norm, ||A||₁ is stored in the upper right corner on output This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error. This function is available only in the C library starting from v2.6.

ldlt-solve(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Vector LDLT.matrix.size1 --> Math::Libgsl::Vector)

This function solves the system Ax = b using the LDT(L) decomposition of A held in the matrix $LDLT which must have been previously computed by ldlt-decomp. In case of error a failure object is returned. This function is available only in the C library starting from v2.6.

ldlt-svx(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Vector LDLT.matrix.size1 --> Int)

This function solves the system Ax = b using the LDT(L) decomposition of A held in the matrix x should contain the right-hand side b, which is replaced by the solution on output. In case of error a failure object is returned. This function is available only in the C library starting from v2.6.

ldlt-rcond(Math::Libgsl::Matrix LDLT.matrix.size2 --> Num)

This function estimates the reciprocal condition number (using the 1-norm) of the symmetric nonsingular matrix A, using its LDT(L) decomposition provided in $LDLT. The function returns the reciprocal condition number estimate, defined as 1/(||A + E||₁ · ||(A + E)⁻¹||₁). In case of error a failure object is returned. This function is available only in the C library starting from v2.6.

symmtd-decomp(Math::Libgsl::Matrix A.matrix.size2 --> Math::Libgsl::Vector)

This function factorizes the symmetric square matrix A contain the tridiagonal matrix T. The remaining lower triangular part of the input matrix contains the Householder vectors which, together with the Householder coefficients tau returned as a Math::Libgsl::Vector object, encode the orthogonal matrix Q. The upper triangular part of $A is not referenced. In case of error a failure object is returned.

symmtd-unpack(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size1 - 1 --> List)

This function unpacks the encoded symmetric tridiagonal decomposition (tau) obtained from symmtd-decomp. The function returns a List: the orthogonal matrix diag as a Math::Libgsl::Vector, and the vector of subdiagonal elements $subdiag as a Math::Libgsl::Vector. In case of error a failure object is returned.

symmtd-unpack-T(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function unpacks the diagonal and subdiagonal of the encoded symmetric tridiagonal decomposition (tau) obtained from symmtd-decomp. The function returns a List of two Math::Libgsl::Vector: subdiag. In case of error a failure object is returned.

hessenberg-decomp(Math::Libgsl::Matrix A.matrix.size2 --> Math::Libgsl::Vector)

This function computes the Hessenberg decomposition of the matrix A. The information required to construct the matrix U is stored in the lower triangular portion of A (below the subdiagonal) and the Householder coefficients are returned as a Math::Libgsl::Vector. In case of error a failure object is returned.

hessenberg-unpack(Math::Libgsl::Matrix H.matrix.size2, Math::Libgsl::Vector H.matrix.size1 --> Math::Libgsl::Matrix)

This function constructs the orthogonal matrix U from the information stored in the Hessenberg matrix tau. tau are outputs from hessenberg-decomp. The function returns the Math::Libgsl::Matrix which stores $U. In case of error a failure object is returned.

hessenberg-unpack-accum(Math::Libgsl::Matrix H.matrix.size2, Math::Libgsl::Vector H.matrix.size1 --> Math::Libgsl::Matrix)

This function is similar to hessenberg-unpack, except it accumulates the matrix U into V, so that V′ = VU. The function returns the Math::Libgsl::Matrix which stores $V. In case of error a failure object is returned.

hessenberg-set-zero(Math::Libgsl::Matrix H.matrix.size2 --> Int)

This function sets the lower triangular portion of $H, below the subdiagonal, to zero. It is useful for clearing out the Householder vectors after calling hessenberg-decomp. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

hesstri-decomp(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Matrix B.matrix.size1 == B.matrix.size2 == similarity = True --> List)

This function computes the Hessenberg-Triangular decomposition of the matrix pair (B). On output, H is stored in B. If the $similarity Bool parameter is True (default), then the similarity transformations are returned as a List of Math::Libgsl::Matrix objects. In case of error a failure object is returned.

bidiag-decomp(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function factorizes the M-by-N matrix A. The orthogonal matrices U and V are stored as compressed Householder vectors in the remaining elements of tau_U and $tau_V. In case of error a failure object is returned.

bidiag-unpack(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size2, Math::Libgsl::Vector A.matrix.size2 - 1 --> List)

This function unpacks the bidiagonal decomposition of A, tau_V) into the separate orthogonal matrices U, V and the diagonal vector diag and superdiagonal superdiag. Note that U is stored as a compact M-by-N orthogonal matrix satisfying T(U)U = I for efficiency. The function returns a List of four objects: the Math::Libgsl::Matrix V, and the Math::Libgsl::Vector sdiag. In case of error a failure object is returned.

bidiag-unpack2(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size2, Math::Libgsl::Vector A.matrix.size2 - 1 --> Math::Libgsl::Matrix)

This function unpacks the bidiagonal decomposition of A, tau_V) into the separate orthogonal matrices U, V and the diagonal vector diag and superdiagonal superdiag. The matrix U is stored in-place in V. In case of error a failure object is returned.

bidiag-unpack-B(Math::Libgsl::Matrix A.matrix.size2 --> List)

This function unpacks the diagonal and superdiagonal of the bidiagonal decomposition of diag and $sdiag. In case of error a failure object is returned.

givens(Num() b --> List)

This function computes c = cos θ and s = sin θ so that the Givens matrix G(θ) acting on the vector (a, b) produces (r, 0), with r = √ a² + b². The function returns a List of Num: the c and s elements of the Givens matrix.

givens-gv(Math::Libgsl::Vector i, Int c, Num() $s)

This function applies the Givens rotation defined by c = cos θ and s = sin θ to the i and j elements of v. On output, (v(i), v(j)) ← G(θ)(v(i), v(j)). This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

householder-transform(Math::Libgsl::Vector $w --> Num)

This function prepares a Householder transformation H = I − τvT(v) which can be used to zero all the elements of the input vector w and the scalar τ is returned. The householder vector v is normalized so that v[0] = 1, however this 1 is not stored in the output vector. Instead, τ.

householder-hm(Num() v, Math::Libgsl::Matrix $A --> Int)

This function applies the Householder matrix H defined by the scalar v to the left-hand side of the matrix A. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

householder-mh(Num() v, Math::Libgsl::Matrix $A --> Int)

This function applies the Householder matrix H defined by the scalar v to the right-hand side of the matrix A. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

householder-hv(Num() v, Math::Libgsl::Vector $w --> Int)

This function applies the Householder transformation H defined by the scalar v to the vector $w. On output the result Hw is stored in w. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

HH-solve(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size1 --> Math::Libgsl::Vector)

This function solves the system Ax = b directly using Householder transformations. A is destroyed by the Householder transformations. The function returns a Math::Libgsl::Vector object: the solution $x. In case of error a failure object is returned.

HH-svx(Math::Libgsl::Matrix A.matrix.size2, Math::Libgsl::Vector A.matrix.size2 --> Int)

This function solves the system Ax = b in-place using Householder transformations. On input A is destroyed by the Householder transformations. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tridiag-solve(Math::Libgsl::Vector abovediag where _.vector.size == belowdiag where _.vector.size == rhs where *.vector.size == $diag.vector.size --> Math::Libgsl::Vector)

This function solves the general N-by-N system Ax = b where A is tridiagonal (N ≥ 2). The super-diagonal and sub-diagonal vectors must be one element shorter than the diagonal vector diag. The function returns the Math::Libgsl::Vector solution $x. In case of error a failure object is returned.

tridiag-symm-solve(Math::Libgsl::Vector offdiag where _.vector.size == rhs where _.vector.size == $diag.vector.size --> Math::Libgsl::Vector)

This function solves the general N-by-N system Ax = b where A is symmetric tridiagonal (N ≥ 2). The off-diagonal vector diag. The function returns the Math::Libgsl::Vector solution $x. In case of error a failure object is returned.

tridiag-cyc-solve(Math::Libgsl::Vector abovediag where _.vector.size == belowdiag where _.vector.size == rhs where _.vector.size == $diag.vector.size --> Math::Libgsl::Vector)

This function solves the general N-by-N system Ax = b where A is cyclic tridiagonal (N ≥ 3). The cyclic super-diagonal and sub-diagonal vectors belowdiag must have the same number of elements as the diagonal vector x. In case of error a failure object is returned.

tridiag-symm-cyc-solve(Math::Libgsl::Vector offdiag where _.vector.size == rhs where *.vector.size == $diag.vector.size --> Math::Libgsl::Vector)

This function solves the general N-by-N system Ax = b where A is symmetric cyclic tridiagonal (N ≥ 3). The cyclic off-diagonal vector diag. The function returns the Math::Libgsl::Vector solution $x. In case of error a failure object is returned.

tri-upper-rcond(Math::Libgsl::Matrix A.matrix.size2 --> Num)

This function estimates the reciprocal condition number. In case of error a failure object is returned.

tri-lower-rcond(Math::Libgsl::Matrix A.matrix.size2 --> Num)

This function estimates the reciprocal condition number. In case of error a failure object is returned.

tri-upper-invert(Math::Libgsl::Matrix T.matrix.size2 --> Int)

This function computes the in-place inverse of the triangular matrix $T, stored in the upper triangle. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-lower-invert(Math::Libgsl::Matrix T.matrix.size2 --> Int)

This function computes the in-place inverse of the triangular matrix $T, stored in the lower triangle. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-upper-unit-invert(Math::Libgsl::Matrix T.matrix.size2 --> Int)

This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-lower-unit-invert(Math::Libgsl::Matrix T.matrix.size2 --> Int)

This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-invert(Int Diag, Math::Libgsl::Matrix T.matrix.size2 --> Int)

This function is available only from the C library v2.6. This function computes the in-place inverse of the triangular matrix Uplo = CblasLower and upper triangle when Diag = CblasUnit, CblasNonUnit specifies whether the matrix is unit triangular. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-LTL(Math::Libgsl::Matrix L.matrix.size2 --> Int)

This function is available only from the C library v2.6. This function computes the product LT(L) in-place and stores it in the lower triangle of $L on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-UL(Math::Libgsl::Matrix LU.matrix.size2 --> Int)

This function is available only from the C library v2.6. This function compute the product LU, as computed by LU-decomp. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-rcond(Int A where *.matrix.size1 == $A.matrix.size2 --> Num)

This function is available only from the C library v2.6. This function estimates the 1-norm reciprocal condition number of the triangular matrix Uplo is CblasLower and upper triangle when $Uplo is CblasUpper. The function returns the reciprocal condition number 1/(||A||₁||A⁻¹||₁). In case of error a failure object is returned.

cholesky-band-decomp(Math::Libgsl::Matrix A.matrix.size2 --> Int)

This function is available only from the C library v2.6. This function factorizes the symmetric, positive-definite square matrix A is given in symmetric banded format, and has dimensions N-by-(p + 1), where p is the lower bandwidth of the matrix. On output, the entries of A is used to store the matrix 1-norm, used later by cholesky-band-rcond() to calculate the reciprocal condition number. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-band-solve(Math::Libgsl::Matrix LLT.matrix.size2, Math::Libgsl::Vector LLT.matrix.size1 --> Math::Libgsl::Vector)

This function is available only from the C library v2.6. This function solves the symmetric banded system Ax = b using the Cholesky decomposition of A held in the matrix x. In case of error a failure object is returned.

cholesky-band-svx(Math::Libgsl::Matrix LLT.matrix.size2, Math::Libgsl::Vector LLT.matrix.size1 --> Int)

This function is available only from the C library v2.6. This function solves the symmetric banded system Ax = b using the Cholesky decomposition of A held in the matrix x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-band-invert(Math::Libgsl::Matrix LLT.matrix.size2 --> Math::Libgsl::Matrix)

This function is available only from the C library v2.6. This function computes the inverse of a symmetric banded matrix from its Cholesky decomposition $LLT, which must have been previously computed by cholesky-band-decomp. The function returns the inverse matrix as a Math::Libgsl::Matrix solution object. In case of error a failure object is returned.

cholesky-band-unpack(Math::Libgsl::Matrix LLT.matrix.size2 --> Math::Libgsl::Matrix)

This function is available only from the C library v2.6. This function unpacks the lower triangular Cholesky factor from $LLT, which returns as a Math::Libgsl::Matrix object. In case of error a failure object is returned.

cholesky-band-rcond(Math::Libgsl::Matrix LLT.matrix.size2 --> Num)

This function is available only from the C library v2.6. This function estimates the reciprocal condition number (using the 1-norm) of the symmetric banded positive definite matrix A, using its Cholesky decomposition provided in $LLT. The reciprocal condition number estimate, defined as 1/(||A||₁ · ||A⁻¹||₁), is returned.

ldlt-band-decomp(Math::Libgsl::Matrix A.matrix.size2 --> Int)

This function is available only from the C library v2.6. This function factorizes the symmetric, non-singular square matrix A is given in symmetric banded format, and has dimensions N-by-(p + 1), where p is the lower bandwidth of the matrix. On output, the entries of $A are replaced by the entries of the matrices D and L in the same format. If the matrix is singular then the decomposition will fail, returning the error code GSL_EDOM . This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

ldlt-band-solve(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Vector LDLT.matrix.size1 --> Math::Libgsl::Vector)

This function is available only from the C library v2.6. This function solves the symmetric banded system Ax = b using the LDT(L) decomposition of A held in the matrix x. In case of error a failure object is returned.

ldlt-band-svx(Math::Libgsl::Matrix LDLT.matrix.size2, Math::Libgsl::Vector LDLT.matrix.size1 --> Int)

This function is available only from the C library v2.6. This function solves the symmetric banded system Ax = b using the LDT(L) decomposition of A held in the matrix x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

ldlt-band-unpack(Math::Libgsl::Matrix LDLT.matrix.size2 --> List)

This function is available only from the C library v2.6. This function unpacks the unit lower triangular factor L from $LDLT. This function returns a Math::Libgsl::Matrix object, which holds the lower triangular portion of the matrix L, and a Math::Libgsl::Vector object, which contains the diagonal matrix D. In case of error a failure object is returned.

ldlt-band-rcond(Math::Libgsl::Matrix LDLT.matrix.size2 --> Num)

This function is available only from the C library v2.6. This function estimates the reciprocal condition number (using the 1-norm) of the symmetric banded nonsingular matrix A, using its LDT(L) decomposition provided in $LDLT. The reciprocal condition number estimate, defined as 1/(||A||₁ · ||A⁻¹||₁), is returned.

balance-matrix(Math::Libgsl::Matrix $A --> Math::Libgsl::Vector)

This function replaces the matrix $A with its balanced counterpart and returns the diagonal elements of the similarity transformation as a Math::Libgsl::Vector object. In case of error a failure object is returned.

Complex

LU-cdecomp(Math::Libgsl::Matrix::Complex64 A.matrix.size2 --> List)

This function factorizes the matrix A into the LU decomposition PA = LU. The factorization is done in place, so on output the diagonal and upper triangular (or trapezoidal) part of the input matrix A contain the matrix U. The lower triangular (or trapezoidal) part of the input matrix (excluding the diagonal) contains L. The diagonal elements of L are unity, and are not stored. The return value is a List: the sign of the permutation and a permutation object, which encodes the permutation matrix P. In case of error a failure object is returned.

LU-csolve(Math::Libgsl::Matrix::Complex64 LU.matrix.size2, Math::Libgsl::Permutation LU.matrix.size1, Math::Libgsl::Vector::Complex64 LU.matrix.size1 --> Math::Libgsl::Vector::Complex64)

This function solves the square system Ax = b using the LU decomposition of A into (LU, p) given by the output of LU-decomp. In case of error a failure object is returned.

LU-csvx(Math::Libgsl::Matrix::Complex64 LU.matrix.size2, Math::Libgsl::Permutation LU.matrix.size1, Math::Libgsl::Vector::Complex64 LU.matrix.size1 --> Int)

This function solves the square system Ax = b in-place using the precomputed LU decomposition of A into (LU, p). On input $x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-crefine(Math::Libgsl::Matrix::Complex64 A.matrix.size2, Math::Libgsl::Matrix::Complex64 LU.matrix.size1 == A.matrix.size1 == p where _.size == b where _.vector.size == x where _.vector.size == $LU.matrix.size1 --> Int)

This function applies an iterative improvement to x, the solution of Ax = b, from the precomputed LU decomposition of A into (LU, p). This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-cinvert(Math::Libgsl::Matrix::Complex64 p where *.size == $LU.matrix.size1 --> Math::Libgsl::Matrix::Complex64)

This function computes the inverse of a matrix A from its LU decomposition (LU, p), returning the matrix inverse. In case of error a failure object is returned.

LU-cdet(Math::Libgsl::Matrix::Complex64 signum where * ~~ -1|1 --> Complex)

This function computes the determinant of a matrix A from its LU decomposition, signum. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-clndet(Math::Libgsl::Matrix::Complex64 $LU --> Num)

This function computes the determinant the logarithm of the absolute value of the determinant of a matrix A, ln |det(A)|, from its LU decomposition, $LU. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

LU-csgndet(Math::Libgsl::Matrix::Complex64 signum where * ~~ -1|1 --> Complex)

This function computes the sign or phase factor of the determinant of a matrix A, det(A)/|det(A)| from its LU decomposition, $LU. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-cdecomp(Math::Libgsl::Matrix::Complex64 A.matrix.size2 --> Int)

This function factorizes the symmetric, positive-definite square matrix A are used (the upper triangular part is ignored). On output the diagonal and lower triangular part of the input matrix $A contain the matrix L, while the upper triangular part contains the original matrix. If the matrix is not positive-definite then the decomposition will fail, returning the error code GSL_EDOM. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-csolve(Math::Libgsl::Matrix::Complex64 A.matrix.size2, Math::Libgsl::Vector::Complex64 A.matrix.size1 --> Math::Libgsl::Vector::Complex64)

These functions solve the system Ax = b using the Cholesky decomposition of x. In case of error a failure object is returned.

cholesky-csvx(Math::Libgsl::Matrix::Complex64 A.matrix.size2, Math::Libgsl::Vector::Complex64 A.matrix.size1 --> Int)

This function solves the system Ax = b in-place using the Cholesky decomposition of x should contain the right-hand side b, which is replaced by the solution on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

cholesky-cinvert(Math::Libgsl::Matrix::Complex64 A.matrix.size2 --> Int)

These functions compute the inverse of a matrix from its Cholesky decomposition A. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

hermtd-cdecomp(Math::Libgsl::Matrix::Complex64 A.matrix.size2 --> Math::Libgsl::Vector::Complex64)

This function factorizes the hermitian matrix A contain the tridiagonal matrix T. The remaining lower triangular part of the input matrix contains the Householder vectors which, together with the Householder coefficients tau returned as a Math::Libgsl::Vector::Complex64 object, encode the orthogonal matrix U. The upper triangular part of $A and the imaginary parts of the diagonal are not referenced. In case of error a failure object is returned.

hermtd-cunpack(Math::Libgsl::Matrix::Complex64 A.matrix.size2, Math::Libgsl::Vector::Complex64 A.matrix.size1 - 1 --> List)

This function unpacks the encoded symmetric tridiagonal decomposition (tau) obtained from hermtd-cdecomp. The function returns a List: the unitary matrix diag as a Math::Libgsl::Vector, and the real vector of subdiagonal elements $subdiag as a Math::Libgsl::Vector. In case of error a failure object is returned.

hermtd-cunpack-T(Math::Libgsl::Matrix::Complex64 A.matrix.size2 --> List)

This function unpacks the diagonal and subdiagonal of the encoded symmetric tridiagonal decomposition (tau) obtained from hermtd-cdecomp. The function returns a List of two real Math::Libgsl::Vector: subdiag. In case of error a failure object is returned.

householder-ctransform(Math::Libgsl::Vector::Complex64 $w --> Complex)

This function prepares a Householder transformation H = I − τvT(v) which can be used to zero all the elements of the input vector w and the scalar τ is returned. The householder vector v is normalized so that v[0] = 1, however this 1 is not stored in the output vector. Instead, τ.

householder-chm(Complex v, Math::Libgsl::Matrix::Complex64 $A --> Int)

This function applies the Householder matrix H defined by the scalar v to the left-hand side of the matrix A. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

householder-cmh(Complex v, Math::Libgsl::Matrix::Complex64 $A --> Int)

This function applies the Householder matrix H defined by the scalar v to the right-hand side of the matrix A. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

householder-chv(Complex v, Math::Libgsl::Vector::Complex64 $w --> Int)

This function applies the Householder transformation H defined by the scalar v to the vector $w. On output the result Hw is stored in w. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-cinvert(Int Diag, Math::Libgsl::Matrix::Complex64 $A --> Int)

This function is available only from the C library v2.6. This function computes the in-place inverse of the triangular matrix Uplo = CblasLower and upper triangle when Diag = CblasUnit, CblasNonUnit specifies whether the matrix is unit triangular. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-cLHL(Math::Libgsl::Matrix::Complex64 $L --> Int)

This function is available only from the C library v2.6. This function computes the product LT(L) in-place and stores it in the lower triangle of $L on output. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

tri-cUL(Math::Libgsl::Matrix::Complex64 $LU --> Int)

This function is available only from the C library v2.6. This function compute the product LU, as computed by LU-cdecomp. This function returns GSL_SUCCESS if successful, or one of the error codes listed in Math::Libgsl::Constants::gsl-error.

C Library Documentation

For more details on libgsl see https://www.gnu.org/software/gsl/. The excellent C Library manual is available here https://www.gnu.org/software/gsl/doc/html/index.html, or here https://www.gnu.org/software/gsl/doc/latex/gsl-ref.pdf in PDF format.

Prerequisites

This module requires the libgsl library to be installed. Please follow the instructions below based on your platform:

Debian Linux and Ubuntu 20.04+

sudo apt install libgsl23 libgsl-dev libgslcblas0

That command will install libgslcblas0 as well, since it's used by the GSL.

Ubuntu 18.04

libgsl23 and libgslcblas0 have a missing symbol on Ubuntu 18.04. I solved the issue installing the Debian Buster version of those three libraries:

Installation

To install it using zef (a module management tool):

$ zef install Math::Libgsl::LinearAlgebra

AUTHOR

Fernando Santagata [email protected]

COPYRIGHT AND LICENSE

Copyright 2020 Fernando Santagata

This library is free software; you can redistribute it and/or modify it under the Artistic License 2.0.

Math::Libgsl::LinearAlgebra v0.0.3

An interface to libgsl, the Gnu Scientific Library - Linear Algebra.

Authors

  • Fernando Santagata

License

Artistic-2.0

Dependencies

Distribution::Builder::MakeFromJSONMath::Libgsl::ConstantsMath::Libgsl::MatrixMath::Libgsl::PermutationMath::Libgsl::BLASMath::Libgsl::Random

Provides

  • Math::Libgsl::LinearAlgebra
  • Math::Libgsl::LinearAlgebra::Complex64
  • Math::Libgsl::Raw::LinearAlgebra

Documentation

The Camelia image is copyright 2009 by Larry Wall. "Raku" is trademark of the Yet Another Society. All rights reserved.