Copyright | David Johnson (c) 2019-2020 |
---|---|
License | BSD 3 |
Maintainer | David Johnson <djohnson.m@gmail.com> |
Stability | Experimental |
Portability | GHC |
Safe Haskell | None |
Language | Haskell2010 |
LAPACK — Linear Algebra PACKage
>>> (u,e,d) = svd (constant @Double [3,3] 10) >>> u ArrayFire Array [3 3 1 1] -0.5774 0.8165 -0.0000 -0.5774 -0.4082 -0.7071 -0.5774 -0.4082 0.7071 >>> e ArrayFire Array [3 1 1 1] 30.0000 0.0000 0.0000 >>> d ArrayFire Array [3 3 1 1] -0.5774 -0.5774 -0.5774 -0.8165 0.4082 0.4082 -0.0000 0.7071 -0.7071
Synopsis
- svd :: AFType a => Array a -> (Array a, Array a, Array a)
- svdInPlace :: AFType a => Array a -> (Array a, Array a, Array a)
- lu :: AFType a => Array a -> (Array a, Array a, Array a)
- luInPlace :: AFType a => Array a -> Bool -> Array a
- qr :: AFType a => Array a -> (Array a, Array a, Array a)
- qrInPlace :: AFType a => Array a -> Array a
- cholesky :: AFType a => Array a -> Bool -> (Int, Array a)
- choleskyInplace :: AFType a => Array a -> Bool -> Int
- solve :: AFType a => Array a -> Array a -> MatProp -> Array a
- solveLU :: AFType a => Array a -> Array a -> Array a -> MatProp -> Array a
- inverse :: AFType a => Array a -> MatProp -> Array a
- rank :: AFType a => Array a -> Double -> Int
- det :: AFType a => Array a -> (Double, Double)
- norm :: AFType a => Array a -> NormType -> Double -> Double -> Double
- isLAPACKAvailable :: Bool
Documentation
:: AFType a | |
=> Array a | the input Matrix |
-> (Array a, Array a, Array a) | Output |
Singular Value Decomposition
The arrayfire function only returns the non zero diagonal elements of S.
:: AFType a | |
=> Array a | the input matrix |
-> (Array a, Array a, Array a) | Output |
Singular Value Decomposition (in-place)
The arrayfire function only returns the non zero diagonal elements of S.
:: AFType a | |
=> Array a | input |
-> Bool | a boolean determining if out is upper or lower triangular |
-> (Int, Array a) | contains the triangular matrix. Multiply |
Perform Cholesky Decomposition
This function decomposes a positive definite matrix A into two triangular matrices.
:: AFType a | |
=> Array a | is the coefficient matrix |
-> Array a | is the measured values |
-> MatProp | determining various properties of matrix a |
-> Array a | is the matrix of unknown variables |
Solve a system of equations
:: AFType a | |
=> Array a | is the output matrix from packed LU decomposition of the coefficient matrix |
-> Array a | is the pivot array from packed LU decomposition of the coefficient matrix |
-> Array a | is the matrix of measured values |
-> MatProp | determining various properties of matrix a |
-> Array a | will contain the matrix of unknown variables |
Solve a system of equations.
:: AFType a | |
=> Array a | is input matrix |
-> Double | is the tolerance value |
-> Int | will contain the rank of in |
Pseudo-inverse
Not implemented in 3.6.4
pinverse :: AFType a => Array a -> Double -> MatProp -> Array a pinverse a d m = op1 a (x y -> af_pinverse x y d (toMatProp m))
Find the rank of the input matrix
This function uses af::qr to find the rank of the input matrix within the given tolerance.
:: AFType a | |
=> Array a | is the input matrix |
-> NormType | specifies the |
-> Double | specifies the value of P when type is one of AF_NORM_VECTOR_P, AF_NORM_MATRIX_L_PQ is used. It is ignored for other values of type |
-> Double | specifies the value of Q when type is AF_NORM_MATRIX_L_PQ. This parameter is ignored if type is anything else |
-> Double | will contain the norm of in |
Find the norm of the input matrix.
This function can return the norm using various metrics based on the type paramter.