| Portability | uses ffi | 
|---|---|
| Stability | provisional | 
| Maintainer | Alberto Ruiz (aruiz at um dot es) | 
Numeric.LinearAlgebra.Algorithms
Contents
Description
High level generic interface to common matrix computations.
Specific functions for particular base types can also be explicitly imported from Numeric.LinearAlgebra.LAPACK.
- class (Product t, Convert t, Container Vector t, Container Matrix t, Normed Matrix t, Normed Vector t) => Field t
- linearSolve :: Field t => Matrix t -> Matrix t -> Matrix t
- luSolve :: Field t => (Matrix t, [Int]) -> Matrix t -> Matrix t
- cholSolve :: Field t => Matrix t -> Matrix t -> Matrix t
- linearSolveLS :: Field t => Matrix t -> Matrix t -> Matrix t
- linearSolveSVD :: Field t => Matrix t -> Matrix t -> Matrix t
- inv :: Field t => Matrix t -> Matrix t
- pinv :: Field t => Matrix t -> Matrix t
- det :: Field t => Matrix t -> t
- invlndet :: (Floating t, Field t) => Matrix t -> (Matrix t, (t, t))
- rank :: Field t => Matrix t -> Int
- rcond :: Field t => Matrix t -> Double
- svd :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)
- fullSVD :: Field t => Matrix t -> (Matrix t, Matrix Double, Matrix t)
- thinSVD :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)
- compactSVD :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)
- singularValues :: Field t => Matrix t -> Vector Double
- leftSV :: Field t => Matrix t -> (Matrix t, Vector Double)
- rightSV :: Field t => Matrix t -> (Vector Double, Matrix t)
- eig :: Field t => Matrix t -> (Vector (Complex Double), Matrix (Complex Double))
- eigSH :: Field t => Matrix t -> (Vector Double, Matrix t)
- eigSH' :: Field t => Matrix t -> (Vector Double, Matrix t)
- eigenvalues :: Field t => Matrix t -> Vector (Complex Double)
- eigenvaluesSH :: Field t => Matrix t -> Vector Double
- eigenvaluesSH' :: Field t => Matrix t -> Vector Double
- geigSH' :: Field t => Matrix t -> Matrix t -> (Vector Double, Matrix t)
- qr :: Field t => Matrix t -> (Matrix t, Matrix t)
- rq :: Field t => Matrix t -> (Matrix t, Matrix t)
- chol :: Field t => Matrix t -> Matrix t
- cholSH :: Field t => Matrix t -> Matrix t
- mbCholSH :: Field t => Matrix t -> Maybe (Matrix t)
- hess :: Field t => Matrix t -> (Matrix t, Matrix t)
- schur :: Field t => Matrix t -> (Matrix t, Matrix t)
- lu :: Field t => Matrix t -> (Matrix t, Matrix t, Matrix t, t)
- luPacked :: Field t => Matrix t -> (Matrix t, [Int])
- expm :: Field t => Matrix t -> Matrix t
- sqrtm :: Field t => Matrix t -> Matrix t
- matFunc :: (Complex Double -> Complex Double) -> Matrix (Complex Double) -> Matrix (Complex Double)
- nullspacePrec :: Field t => Double -> Matrix t -> [Vector t]
- nullVector :: Field t => Matrix t -> Vector t
- nullspaceSVD :: Field t => Either Double Int -> Matrix t -> (Vector Double, Matrix t) -> [Vector t]
- class RealFloat (RealOf t) => Normed c t where
- data NormType
- relativeError :: (Normed c t, Container c t) => c t -> c t -> Int
- eps :: Double
- peps :: RealFloat x => x
- i :: Complex Double
- haussholder :: Field a => a -> Vector a -> Matrix a
- unpackQR :: Field t => (Matrix t, Vector t) -> (Matrix t, Matrix t)
- unpackHess :: Field t => (Matrix t -> (Matrix t, Vector t)) -> Matrix t -> (Matrix t, Matrix t)
- pinvTol :: Double -> Matrix Double -> Matrix Double
- ranksv :: Double -> Int -> [Double] -> Int
- full :: (Num t1, Storable t1) => (Matrix t -> (t2, Vector t1, t3)) -> Matrix t -> (t2, Matrix t1, t3)
- economy :: (Element t, Element t1, Element t2) => (Matrix t -> (Matrix t1, Vector Double, Matrix t2)) -> Matrix t -> (Matrix t1, Vector Double, Matrix t2)
Supported types
class (Product t, Convert t, Container Vector t, Container Matrix t, Normed Matrix t, Normed Vector t) => Field t Source
Linear Systems
linearSolve :: Field t => Matrix t -> Matrix t -> Matrix tSource
Solve a linear system (for square coefficient matrix and several right-hand sides) using the LU decomposition. For underconstrained or overconstrained systems use linearSolveLS or linearSolveSVD.
 It is similar to luSolve . luPacked, but linearSolve raises an error if called on a singular system.
luSolve :: Field t => (Matrix t, [Int]) -> Matrix t -> Matrix tSource
Solution of a linear system (for several right hand sides) from the precomputed LU factorization obtained by luPacked.
cholSolve :: Field t => Matrix t -> Matrix t -> Matrix tSource
Solve a symmetric or Hermitian positive definite linear system using a precomputed Cholesky decomposition obtained by chol.
linearSolveLS :: Field t => Matrix t -> Matrix t -> Matrix tSource
Least squared error solution of an overconstrained linear system, or the minimum norm solution of an underconstrained system. For rank-deficient systems use linearSolveSVD.
linearSolveSVD :: Field t => Matrix t -> Matrix t -> Matrix tSource
Minimum norm solution of a general linear least squares problem Ax=B using the SVD. Admits rank-deficient systems but it is slower than linearSolveLS. The effective rank of A is determined by treating as zero those singular valures which are less than eps times the largest singular value.
det :: Field t => Matrix t -> tSource
Determinant of a square matrix. To avoid possible overflow or underflow use invlndet.
Arguments
| :: (Floating t, Field t) | |
| => Matrix t | |
| -> (Matrix t, (t, t)) | (inverse, (log abs det, sign or phase of det)) | 
Joint computation of inverse and logarithm of determinant of a square matrix.
rcond :: Field t => Matrix t -> DoubleSource
Reciprocal of the 2-norm condition number of a matrix, computed from the singular values.
Matrix factorizations
Singular value decomposition
svd :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)Source
Full singular value decomposition.
fullSVD :: Field t => Matrix t -> (Matrix t, Matrix Double, Matrix t)Source
A version of svd which returns an appropriate diagonal matrix with the singular values.
If (u,d,v) = fullSVD m then m == u <> d <> trans v.
thinSVD :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)Source
A version of svd which returns only the min (rows m) (cols m) singular vectors of m.
If (u,s,v) = thinSVD m then m == u <> diag s <> trans v.
compactSVD :: Field t => Matrix t -> (Matrix t, Vector Double, Matrix t)Source
Similar to thinSVD, returning only the nonzero singular values and the corresponding singular vectors.
leftSV :: Field t => Matrix t -> (Matrix t, Vector Double)Source
Singular values and all right singular vectors.
rightSV :: Field t => Matrix t -> (Vector Double, Matrix t)Source
Singular values and all right singular vectors.
Eigensystems
eig :: Field t => Matrix t -> (Vector (Complex Double), Matrix (Complex Double))Source
Eigenvalues and eigenvectors of a general square matrix.
If (s,v) = eig m then m <> v == v <> diag s
eigSH :: Field t => Matrix t -> (Vector Double, Matrix t)Source
Eigenvalues and Eigenvectors of a complex hermitian or real symmetric matrix.
If (s,v) = eigSH m then m == v <> diag s <> ctrans v
eigSH' :: Field t => Matrix t -> (Vector Double, Matrix t)Source
Similar to eigSH without checking that the input matrix is hermitian or symmetric. It works with the upper triangular part.
eigenvalues :: Field t => Matrix t -> Vector (Complex Double)Source
Eigenvalues of a general square matrix.
eigenvaluesSH :: Field t => Matrix t -> Vector DoubleSource
Eigenvalues of a complex hermitian or real symmetric matrix.
eigenvaluesSH' :: Field t => Matrix t -> Vector DoubleSource
Similar to eigenvaluesSH without checking that the input matrix is hermitian or symmetric. It works with the upper triangular part.
Generalized symmetric positive definite eigensystem Av = lBv, for A and B symmetric, B positive definite (conditions not checked).
QR
qr :: Field t => Matrix t -> (Matrix t, Matrix t)Source
QR factorization.
If (q,r) = qr m then m == q <> r, where q is unitary and r is upper triangular.
rq :: Field t => Matrix t -> (Matrix t, Matrix t)Source
RQ factorization.
If (r,q) = rq m then m == r <> q, where q is unitary and r is upper triangular.
Cholesky
chol :: Field t => Matrix t -> Matrix tSource
Cholesky factorization of a positive definite hermitian or symmetric matrix.
If c = chol m then c is upper triangular and m == ctrans c <> c.
cholSH :: Field t => Matrix t -> Matrix tSource
Similar to chol, without checking that the input matrix is hermitian or symmetric. It works with the upper triangular part.
Hessenberg
hess :: Field t => Matrix t -> (Matrix t, Matrix t)Source
Hessenberg factorization.
If (p,h) = hess m then m == p <> h <> ctrans p, where p is unitary
 and h is in upper Hessenberg form (it has zero entries below the first subdiagonal).
Schur
schur :: Field t => Matrix t -> (Matrix t, Matrix t)Source
Schur factorization.
If (u,s) = schur m then m == u <> s <> ctrans u, where u is unitary
 and s is a Shur matrix. A complex Schur matrix is upper triangular. A real Schur matrix is
 upper triangular in 2x2 blocks.
"Anything that the Jordan decomposition can do, the Schur decomposition can do better!" (Van Loan)
LU
lu :: Field t => Matrix t -> (Matrix t, Matrix t, Matrix t, t)Source
Explicit LU factorization of a general matrix.
If (l,u,p,s) = lu m then m == p <> l <> u, where l is lower triangular,
 u is upper triangular, p is a permutation matrix and s is the signature of the permutation.
luPacked :: Field t => Matrix t -> (Matrix t, [Int])Source
Obtains the LU decomposition of a matrix in a compact data structure suitable for luSolve.
Matrix functions
expm :: Field t => Matrix t -> Matrix tSource
Matrix exponential. It uses a direct translation of Algorithm 11.3.1 in Golub & Van Loan, based on a scaled Pade approximation.
sqrtm :: Field t => Matrix t -> Matrix tSource
Matrix square root. Currently it uses a simple iterative algorithm described in Wikipedia.
It only works with invertible matrices that have a real solution. For diagonalizable matrices you can try matFunc sqrt.
m = (2><2) [4,9
           ,0,4] :: Matrix Double>sqrtm m (2><2) [ 2.0, 2.25 , 0.0, 2.0 ]
matFunc :: (Complex Double -> Complex Double) -> Matrix (Complex Double) -> Matrix (Complex Double)Source
Generic matrix functions for diagonalizable matrices. For instance:
logm = matFunc log
Nullspace
Arguments
| :: Field t | |
| => Double | |
| -> Matrix t | input matrix | 
| -> [Vector t] | list of unitary vectors spanning the nullspace | 
The nullspace of a matrix. See also nullspaceSVD.
nullVector :: Field t => Matrix t -> Vector tSource
The nullspace of a matrix, assumed to be one-dimensional, with machine precision.
Arguments
| :: Field t | |
| => Either Double Int | Left "numeric" zero (eg. 1* | 
| -> Matrix t | input matrix m | 
| -> (Vector Double, Matrix t) | 
 | 
| -> [Vector t] | list of unitary vectors spanning the nullspace | 
The nullspace of a matrix from its SVD decomposition.
Norms
relativeError :: (Normed c t, Container c t) => c t -> c t -> IntSource
Approximate number of common digits in the maximum element.
Misc
The machine precision of a Double: eps = 2.22044604925031e-16 (the value used by GNU-Octave).
Util
haussholder :: Field a => a -> Vector a -> Matrix aSource
unpackHess :: Field t => (Matrix t -> (Matrix t, Vector t)) -> Matrix t -> (Matrix t, Matrix t)Source
Arguments
| :: Double | numeric zero (e.g. 1* | 
| -> Int | maximum dimension of the matrix | 
| -> [Double] | singular values | 
| -> Int | rank of m | 
Numeric rank of a matrix from its singular values.