The S-lang interpreter has wide-spread native support for
manipulating arrays and matrices. The gslmatrix
supplements
this by adding some linear algebra routines such LU decomposition as
well as routines for dealing with eigenvalues and eigenvectors.
GSL has separate functions for complex numbers. Rather than creating separater wrappers for each of these functions, the complex-valued routines have been incorporated into single wrappers that support for both real and complex numbers. In this way the interface is polymorphic.
Factorize a square matrix into its LU decomposition
(LU,p) = linalg_LU_decomp (A [,&signum])
This routines returns the LU decomposition of the square matrix
A
such that P#A == LU
. See the corresponding GSL
documentation for how L
and U
are stored in LU
,
and how the permutation matrix P
is defined. For many
applications, it is unnecessary to unpack the matrix LU
into
its separate components.
If the optional argument &signum
is given, upon return
signum
will be set to the sign of the permutation that relates
P
to the identity matrix.
Compute the determinant of a matrix from its LU decomposition
det = linalg_LU_det (LU, signum)
This function computes the determinant of a matrix from its LU decomposition. In the LU form, determinant is given by the product of the diagonal elements with the sign of the permutation.
require ("gslmatrix");
define determinant (A)
{
variable LU, sig;
(LU,) = linalg_LU_decomp (A, &sig);
return linalg_LU_det (LU,sig);
}
linalg_LU_lndet,
linalg_LU_decomp,
linalg_LU_invert,
linalg_LU_solve
Compute the log of a determinant using LU decomposition
det = linalg_LU_lndet (LU)
This function computes the natural logarithm of the determinant of a matrix from its LU decomposition. In the LU form, determinant is given by the product of the diagonal elements with the sign of the permutation. This function is useful for cases where the product of the diagonal elements would overflow.
linalg_LU_det,
linalg_LU_decomp,
linalg_LU_solve,
linalg_LU_invert
Compute the inverse of a matrix via its LU decomposition
inv = linalg_LU_invert (LU, p)
This function may be used to compute the inverse of a matrix from
its LU decomposition. For the purposes of inverting a set of linear
equations, it is preferable to use the linalg_LU_solve
function rather than inverting the equations via the inverse.
define matrix_inverse (A)
{
return linalg_LU_invert (linalg_LU_decomp (A));
}
Solve a set of linear equations using LU decomposition
x = linalg_LU_solve (LU, p, b)
This function solves the square linear system of equations
A#x=b
for the vector x
via the LU decomposition of
A
.
define solve_equations (A, b)
{
return linalg_LU_solve (linalg_LU_decomp (A), b);
}
Factor a matrix into its QR form
(QR, tau) = linalg_QR_decomp(A)
This function may be used to decompose a rectangular matrix into its
so-called QR such that A=Q#R
where Q
is a square
orthogonal matrix and R
is a rectangular right-triangular
matrix.
The factor R
encoded in the diagonal and
upper-triangular elements of the first return value QR
. The
matrix Q
is encoded in the lower triangular part of QR
and the vector tau
via Householder vectors and coefficients.
See the corresponding
GSL documentation for the details of the
encoding. For most uses encoding details are not required.
Solve a system of linear equations using QR decomposition
x = linalg_QR_solve(QR, tau, b [,&residual])
This function may be used to solve the linear system A#x=b
using the QR
decomposition of A
.
If the optional fourth argument is present (&residual
), or if
QR
is not a square matrix, then the linear system will be
solved in the least-squares sense by minimizing the (Euclidean) norm
of A#x-b
. Upon return, the value of the variable
residual
is set to the the norm of A#x-b
.
GSL has a separate function called gsl_linalg_QR_lssolve
for
computing this least-squares solution. The linalg_QR_solve
combines both gsl_linalg_QR_lssolve
and
gsl_linalg_QR_solve
into a single routine.
Perform a singular-value decomposition on a matrix
(U,S,V) = linalg_SV_decomp(A)
This function factors a MxN (M>=N) rectangular matrix A
into
three factors such that A = U#S#transpose(V)
, where S
is diagonal matrix containing the singular values of A
and
V
is a square orthogonal matrix. Since S
is diagonal,
it is returned as a 1-d array.
Solve a linear system using Singular-Value Decomposition
x = linalg_SV_solve (U,V,S,b)
This function ``solves'' the linear system A#x=b
using the
SVD form of A
.
define svd_solve (A, b)
{
variable U, V, S;
(U,V,S) = linalg_SV_decomp (A);
return linalg_SV_solve (U,V,S,b);
}
Compute the eigenvalues and eigenvectors of a Hermitian matrix
(eigvecs, eigvals)=eigen_symmv(A)
This function computes the eigenvalues and eigenvectors of a
Hermitian (or real-symmetric) square matrix A
. The
eigenvalues are returned sorted on their absolute value (or norm) in
descending order.
Compute the eigenvalues and eigenvectors of a matrix
(eigvecs, eigvals)=eigen_nonsymmv(A)
This function returns the eigenvalues and eigenvectors of a real
non-symmetric matrix A
. As such quantities are in general
complex, complex-valued arrays will be returned. The eigenvalues
are returned in descending order sorted upon norm.