Title: Linear Algebra

1 Module/Classes

2 LU Decomposition

GSL::Linalg::LU.decomp(A)
GSL::Matrix#LU_decomp
These method calculate the LU decomposition of the matrix. The returned value is an array of [LU, perm, sign].
GSL::Linalg::LU.solve(A, b)
GSL::Linalg::LU.solve(lu, perm, b)
GSL::Matrix#LU_solve(b)
GSL::Linalg::LU::LUMatrix#solve(perm, b)

The following is an example to solve a linear system

A x = b, b = [1, 2, 3, 4]

using LU decomposition.

ex1)

require("gsl")
A = Matrix.new([0.18, 0.60, 0.57, 0.96], [0.41, 0.24, 0.99, 0.58],
                   [0.14, 0.30, 0.97, 0.66], [0.51, 0.13, 0.19, 0.85])
lu, perm, sign = A.LU_decomp         
p lu.class                         <- GSL::Linalg::LU::LUMatrix, subclass of the Matrix class
b = [1, 2, 3, 4].to_gv             <- "b" must be given as a GSL::Vector
x = Linalg::LU.solve(lu, perm, b)    
p x                                <- The solution is returned by a vector.

ex2)

require("gsl")
A = Matrix.new([0.18, 0.60, 0.57, 0.96], [0.41, 0.24, 0.99, 0.58],
                   [0.14, 0.30, 0.97, 0.66], [0.51, 0.13, 0.19, 0.85])
p A.class                           <- GSL::Matrix
x = Linalg::LU.solve(A, b)          <- LU decomposition is calculated internally (A is not modified)
GSL::Linalg::LU.svx(A, b)
GSL::Linalg::LU.svx(lu, perm, b)
GSL::Matrix#svx(b)
GSL::Linalg::LU::LUMatrix#svx(perm, b)
These solve the system Ax = b. The input vector b is overwitten by the solution x.
GSL::Linalg::LU.refine(A, lu, perm, b, x)
This method applys an iterative improvement to x, the solution of A x = b, using the LU decomposition of A.
GSL::Linalg::LU.invert(A)
GSL::Linalg::LU.invert(lu, perm)
GSL::Matrix#invert
GSL::Linalg::LU::LUMatrix#invert(perm)
These computes and returns the inverse of the matrix.
GSL::Linalg::LU.det(A)
GSL::Linalg::LU.det(lu, signum)
GSL::Matrix#det
GSL::Linalg::LU::LUMatrix#det(signum)
These methods return the determinant of the matrix.

3 Complex LU decomposition

GSL::Linalg::Complex::LU_decomp!(A)
GSL::Linalg::Complex::LU::decomp!(A)
GSL::Matrix::Complex#LU_decomp!
GSL::Matrix::Complex#LU_decomp!
Factorizes the square matrix A into the LU decomposition PA = LU, and returns an array, [perm, signum]. A is changed.
GSL::Linalg::Complex::LU_decomp(A)
GSL::Linalg::Complex::LU::decomp(A)
GSL::Matrix::Complex#LU_decomp
Factorizes the square matrix A into the LU decomposition PA = LU, and returns an array, [LU, perm, signum]. A is not changed.
GSL::Linalg::Complex::LU_solve(A, b)
GSL::Linalg::Complex::LU::solve(A, b)
GSL::Linalg::Complex::LU_solve(A, b)
GSL::Matrix::Complex#LU_solve(b)
GSL::Linalg::Complex::solve(LU, perm, b)
GSL::Linalg::Complex::LU::solve(LU, perm, b)
GSL::Linalg::Complex::LU::LUMatirx#solve(perm, b)
GSL::Linalg::Complex::LU_svx(A, x)
GSL::Linalg::Complex::LU::svx(A, x)
GSL::Linalg::Complex::LU_svx(A, x)
GSL::Matrix::Complex#LU_svx(x)
GSL::Linalg::Complex::svx(LU, perm, x)
GSL::Linalg::Complex::LU::svx(LU, perm, x)
GSL::Linalg::Complex::LU::LUMatirx#svx(perm, x)
GSL::Linalg::Complex::LU_refine(A, LU, perm, b, x)
GSL::Linalg::Complex::LU_::refine(A, LU, perm, b, x)
GSL::Linalg::Complex::LU_invert(A)
GSL::Linalg::Complex::LU::invert(A)
GSL::Linalg::Complex::LU_invert(LU, perm)
GSL::Linalg::Complex::LU::invert(LU, perm)
GSL::Matrix::Complex#LU_invert
GSL::Matrix::Complex#invert
GSL::Linalg::Complex::LU::LUMatrix#invert(perm)
GSL::Linalg::Complex::LU_det(A)
GSL::Linalg::Complex::LU::det(A)
GSL::Linalg::Complex::LU_det(LU, signum)
GSL::Linalg::Complex::LU::det(LU, signum)
GSL::Matrix::Complex#LU_det
GSL::Matrix::Complex#det
GSL::Linalg::Complex::LU::LUMatrix#det(signum)
GSL::Linalg::Complex::LU_lndet(A)
GSL::Linalg::Complex::LU::lndet(A)
GSL::Linalg::Complex::LU_lndet(LU)
GSL::Linalg::Complex::LU::lndet(LU)
GSL::Matrix::Complex#LU_lndet
GSL::Matrix::Complex#lndet
GSL::Linalg::Complex::LU::LUMatrix#lndet
GSL::Linalg::Complex::LU_sgndet(A)
GSL::Linalg::Complex::LU::sgndet(A)
GSL::Linalg::Complex::LU_sgndet(LU, signum)
GSL::Linalg::Complex::LU::sgndet(LU, signum)
GSL::Matrix::Complex#LU_sgndet
GSL::Matrix::Complex#sgndet
GSL::Linalg::Complex::LU::LUMatrix#sgndet(signum)

4 QR decomposition

GSL::Linalg::QR.decomp(A)
GSL::Matrix#QR_decomp

These compute QR decomposition of the matrix and return an array [QR, tau].

or

qr, tau = QR.decomp(m)
GSL::Linalg::QR.solve(A, b)
GSL::Linalg::QR.solve(QR, tau, b)
GSL::Matrix#QR_solve(b)
GSL::Linalg::QR::QRMatrix#solve(tau, b)
Solve the system A x = b using the QR decomposition.
GSL::Linalg::QR.svx(A, x)
GSL::Linalg::QR.svx(QR, tau, x)
GSL::Matrix#QR_svx(x)
GSL::Linalg::QR::QRMatrix#svx(tau, x)
Solve the system A x = b. The input vector x is first give by the right-hand side vector b, and is overwritten by the solution.
GSL::Linalg::QR.unpack(QR, tau)
GSL::Linalg::QR::QRMatrix#unpack(tau)
Unpack the encoded QR decomposition QR,tau and return an array [Q, R].
GSL::Linalg::QR.QRsolve(Q, R, tau)
This method solves the system R x = Q^T b for x. It can be used when the QR decomposition of a matrix is available in unpacked form as Q,R.

5 QR Decomposition with Column Pivoting

GSL::Linalg::QRPT.decomp(A)
GSL::Matrix#QRPT_decomp
These methods factorize the M-by-N matrix A into the QRP^T decomposition A = Q R P^T, and return an array [QR, tau, perm, signum].
GSL::Linalg::QRPT.decomp2(A)
GSL::Matrix#QRPT_decomp2
These return an array [Q, R, tau, perm, signum].
GSL::Linalg::QRPT.solve(m, b)
GSL::Linalg::QRPT.solve(qr, tau, perm, b)
GSL::Matrix#QRPT_solve(A, b)
GSL::Linalg::QRPT::QRMatrix#solve(qr, tau, perm, b)
These methods solve the system A x = b using the QRP^T decomposition of A into QR, tau, perm. The solution x is returned as a Vector.
GSL::Linalg::QRPT.svx(m, b)
GSL::Linalg::QRPT.svx(qr, tau, perm, b)
GSL::Matrix#QRPT_svx(A, b)
These methods solve the system A x = b using the QRP^T decomposition of A into QR, tau, perm. The input b is overwritten by the solution x.
GSL::Linalg::QRPT.QRsolve(q, r, tau, perm, b)
This method solves the system R P^T x = Q^T b for x. It can be used when the QR decomposition of a matrix is available in unpacked form as q, r obtained by the method decomp2.
GSL::Linalg::QRPT.update(q, r, perm, u, v)
GSL::Linalg::QRPT.Rsolve(qr, perm, b)
GSL::Linalg::QRPT::QRMatrix#Rsolve(perm, b)
GSL::Linalg::QRPT.Rsvx(qr, perm, b)
GSL::Linalg::QRPT::QRMatrix#Rsvx(perm, b)

6 Singular Value Decomposition

GSL::Linalg::SV.decomp(A)
GSL::Matrix#SV_decomp
These methods factorize the M-by-N matrix A into the singular value decomposition A = U S V^T using the Golub-Reinsch SVD algorithm, and return an array [U, V, S].
GSL::Linalg::SV.decomp_mod(A)
GSL::Matrix#SV_decomp_mod
These compute the SVD using the modified Golub-Reinsch algorithm, which is faster for M>>N.
GSL::Linalg::SV.decomp_jacobi(A)
GSL::Matrix#SV_decomp_jacobi
These compute the SVD using one-sided Jacobi orthogonalization. The Jacobi method can compute singular values to higher relative accuracy than Golub-Reinsch algorithms.
GSL::Linalg::SV.solve(A, b)
GSL::Linalg::SV.solve(U, V, S, b)
GSL::Matrix#SV_solve(b)
These methods solve the system A x = b using the singular value decomposition U, S, V of A.

7 Bidiagonalization

GSL::Linalg::Bidiag::decomp!(A)
GSL::Linalg::bidiag_decomp!(A)
GSL::Linalg::Bidiag::decomp(A)
GSL::Linalg::bidiag_decomp(A)
GSL::Linalg::Bidiag::unpack
GSL::Linalg::bidiag_unpack
GSL::Linalg::Bidiag::unpack2
GSL::Linalg::bidiag_unpack2
GSL::Linalg::Bidiag::unpack_B
GSL::Linalg::bidiag_unpack_B

8 Householder Transformations

GSL::Linalg::Householder::transform(v)
GSL::Linalg::HH::transform(v)
GSL::Vector#householder_transform
These methods prepare a Householder transformation P = I - tau v v^T which can be used to zero all the elements of the input vector except the first. On output the transformation is stored in the vector v and the scalar tau is returned.
GSL::Linalg::Householder::hm(tau, v, A)
GSL::Linalg::HH::hm(tau, v, A)
These methods apply the Householder matrix P defined by the scalar tau and the vector v to the left-hand side of the matrix A. On output the result P A is stored in A.
GSL::Linalg::Householder::mh(tau, v, A)
GSL::Linalg::HH::mh(tau, v, A)
These methods apply the Householder matrix P defined by the scalar tau and the vector v to the right-hand side of the matrix A. On output the result A P is stored in A.
GSL::Linalg::Householder::hv(tau, v, w)
GSL::Linalg::HH::hv(tau, v, w)
These methods apply the Householder transformation P defined by the scalar tau and the vector v to the vector w. On output the result P w is stored in w.

9 Householder solver for linear systems

GSL::Linalg::HH.solve(A, b)
GSL::Matrix#HH_solve(b)
These methods solve the system A x = b directly using Householder transformations.
GSL::Linalg::HH.svx(A, b)
GSL::Matrix#HH_svx(b)
These methods solve the system A x = b in-place directly using Householder transformations. The input vector b is replaced by the solution.

10 Tridiagonal Systems

GSL::Linglg::solve_tridiag(diag, e, f, b)
GSL::Linglg::Tridiag::solve(diag, e, f, b)

These methods solve the general N-by-N system A x = b where A is tridiagonal ( N >= 2). The super-diagonal and sub-diagonal vectors e and f must be one element shorter than the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A = ( d_0 e_0  0   0  )
    ( f_0 d_1 e_1  0  )
    (  0  f_1 d_2 e_2 )
    (  0   0  f_2 d_3 )
GSL::Linglg::solve_symm_tridiag(diag, e, b)
GSL::Linglg::Tridiag::solve_symm(diag, e, b)

These methods solve the general N-by-N system A x = b where A is symmetric tridiagonal ( N >= 2). The off-diagonal vector e must be one element shorter than the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A = ( d_0 e_0  0   0  )
    ( e_0 d_1 e_1  0  )
    (  0  e_1 d_2 e_2 )
    (  0   0  e_2 d_3 )
GSL::Linglg::solve_cyc_tridiag(diag, e, f, b)
GSL::Linglg::Tridiag::solve_cyc(diag, e, f, b)

These methods solve the general N-by-N system A x = b where A is cyclic tridiagonal ( N >= 3). The cyclic super-diagonal and sub-diagonal vectors e and f must have the same number of elements as the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A = ( d_0 e_0  0  f_3 )
    ( f_0 d_1 e_1  0  )
    (  0  f_1 d_2 e_2 )
    ( e_3  0  f_2 d_3 )
GSL::Linglg::solve_symm_cyc_tridiag(diag, e, b)
GSL::Linglg::Tridiag::solve_symm_cyc(diag, e, b)

These methods solve the general N-by-N system A x = b where A is symmetric cyclic tridiagonal ( N >= 3). The cyclic off-diagonal vector e must have the same number of elements as the diagonal vector diag. The form of A for the 4-by-4 case is shown below,

A = ( d_0 e_0  0  e_3 )
    ( e_0 d_1 e_1  0  )
    (  0  e_1 d_2 e_2 )
    ( e_3  0  e_2 d_3 )

11 Balancing

GSL::Linalg.balance_columns(A, D)
GSL::Matrix#balance_columns(D)

back