|
Numeric.GSL.Matrix | Portability | portable (uses FFI) | Stability | provisional | Maintainer | Alberto Ruiz <aruiz@um.es> |
|
|
|
Description |
A few linear algebra computations based on the Numeric.GSL (http://www.gnu.org/software/Numeric.GSL).
|
|
Synopsis |
|
eigSg :: Matrix Double -> (Vector Double, Matrix Double) | | eigHg :: Matrix (Complex Double) -> (Vector Double, Matrix (Complex Double)) | | svdg :: Matrix Double -> (Matrix Double, Vector Double, Matrix Double) | | qr :: Matrix Double -> (Matrix Double, Matrix Double) | | qrPacked :: Matrix Double -> (Matrix Double, Vector Double) | | unpackQR :: (Matrix Double, Vector Double) -> (Matrix Double, Matrix Double) | | cholR :: Matrix Double -> Matrix Double | | cholC :: Matrix (Complex Double) -> Matrix (Complex Double) | | luSolveR :: Matrix Double -> Matrix Double -> Matrix Double | | luSolveC :: Matrix (Complex Double) -> Matrix (Complex Double) -> Matrix (Complex Double) | | luR :: Matrix Double -> (Matrix Double, Matrix Double, [Int], Double) | | luC :: Matrix (Complex Double) -> (Matrix (Complex Double), Matrix (Complex Double), [Int], Complex Double) |
|
|
Documentation |
|
|
eigendecomposition of a real symmetric matrix using gsl_eigen_symmv.
> let (l,v) = eigS $ 'fromLists' [[1,2],[2,1]]
> l
3.000 -1.000
> v
0.707 -0.707
0.707 0.707
> v <> diag l <> trans v
1.000 2.000
2.000 1.000
|
|
|
eigendecomposition of a complex hermitian matrix using gsl_eigen_hermv
> let (l,v) = eigH $ 'fromLists' [[1,2+i],[2-i,3]]
> l
4.449 -0.449
> v
-0.544 0.839
(-0.751,0.375) (-0.487,0.243)
> v <> diag l <> (conjTrans) v
1.000 (2.000,1.000)
(2.000,-1.000) 3.000
|
|
|
Singular value decomposition of a real matrix, using gsl_linalg_SV_decomp_mod:
> let (u,s,v) = svdg $ fromLists [[1,2,3],[-4,1,7]]
\
> u
0.310 -0.951
0.951 0.310
\
> s
8.497 2.792
\
> v
-0.411 -0.785
0.185 -0.570
0.893 -0.243
\
> u <> diag s <> trans v
1. 2. 3.
-4. 1. 7. |
|
|
QR decomposition of a real matrix using gsl_linalg_QR_decomp and gsl_linalg_QR_unpack.
> let (q,r) = qr $ fromLists [[1,3,5,7],[2,0,-2,4]]
\
> q
-0.447 -0.894
-0.894 0.447
\
> r
-2.236 -1.342 -0.447 -6.708
0. -2.683 -5.367 -4.472
\
> q <> r
1.000 3.000 5.000 7.000
2.000 0. -2.000 4.000 |
|
|
|
|
|
|
Cholesky decomposition of a symmetric positive definite real matrix using gsl_linalg_cholesky_decomp.
> chol $ (2><2) [1,2,
2,9::Double]
(2><2)
[ 1.0, 0.0
, 2.0, 2.23606797749979 ] |
|
|
|
|
|
|
|
|
The LU decomposition of a square matrix. Is based on gsl_linalg_LU_decomp and gsl_linalg_complex_LU_decomp as described in http://www.gnu.org/software/Numeric.GSL/manual/Numeric.GSL-ref_13.html#SEC223.
> let m = fromLists [[1,2,-3],[2+3*i,-7,0],[1,-i,2*i]]
> let (l,u,p,s) = luR m L is the lower triangular:
> l
1. 0. 0.
0.154-0.231i 1. 0.
0.154-0.231i 0.624-0.522i 1. U is the upper triangular:
> u
2.+3.i -7. 0.
0. 3.077-1.615i -3.
0. 0. 1.873+0.433i p is a permutation:
> p
[1,0,2] L * U obtains a permuted version of the original matrix:
> extractRows p m
2.+3.i -7. 0.
1. 2. -3.
1. -1.i 2.i
> l <> u
2.+3.i -7. 0.
1. 2. -3.
1. -1.i 2.i s is the sign of the permutation, required to obtain sign of the determinant:
> s * product (toList $ takeDiag u)
(-18.0) :+ (-16.000000000000004)
> LinearAlgebra.Algorithms.det m
(-18.0) :+ (-16.000000000000004) |
|
|
Complex version of luR.
|
|
Produced by Haddock version 2.4.2 |