arrayfire-0.5.0.0: Haskell bindings to the ArrayFire general-purpose GPU library

CopyrightDavid Johnson (c) 2019-2020
LicenseBSD 3
MaintainerDavid Johnson <djohnson.m@gmail.com>
StabilityExperimental
PortabilityGHC
Safe HaskellNone
LanguageHaskell2010

ArrayFire.LAPACK

Description

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

Documentation

svd Source #

Arguments

:: AFType a 
=> Array a

the input Matrix

-> (Array a, Array a, Array a)

Output Array containing (U, diagonal values of sigma, V^H)

Singular Value Decomposition

ArrayFire Docs

The arrayfire function only returns the non zero diagonal elements of S.

svdInPlace Source #

Arguments

:: AFType a 
=> Array a

the input matrix

-> (Array a, Array a, Array a)

Output Array containing (U, diagonal values of sigma, V^H)

Singular Value Decomposition (in-place)

ArrayFire Docs

The arrayfire function only returns the non zero diagonal elements of S.

lu Source #

Arguments

:: AFType a 
=> Array a

is the input matrix

-> (Array a, Array a, Array a)

Returns the output Arrays (lower, upper, pivot)

Perform LU decomposition

ArrayFire Docs

C Interface for LU decomposition.

luInPlace Source #

Arguments

:: AFType a 
=> Array a

contains the input on entry, the packed LU decomposition on exit.

-> Bool

specifies if the pivot is returned in original LAPACK compliant format

-> Array a

will contain the permutation indices to map the input to the decomposition

Perform LU decomposition (in-place).

ArrayFire Docs

C Interface for in place LU decomposition.

qr Source #

Arguments

:: AFType a 
=> Array a

the input matrix

-> (Array a, Array a, Array a)

Returns (q, r, tau) Arrays q is the orthogonal matrix from QR decomposition r is the upper triangular matrix from QR decomposition tau will contain additional information needed for solving a least squares problem using q and r

Perform QR decomposition

ArrayFire Docs

C Interface for QR decomposition.

qrInPlace Source #

Arguments

:: AFType a 
=> Array a

is the input matrix on entry. It contains packed QR decomposition on exit

-> Array a

will contain additional information needed for unpacking the data

Perform QR decomposition

ArrayFire Docs

C Interface for QR decomposition.

cholesky Source #

Arguments

:: AFType a 
=> Array a

input Array

-> Bool

a boolean determining if out is upper or lower triangular

-> (Int, Array a)

contains the triangular matrix. Multiply Int with its conjugate transpose reproduces the input array. is 0 if cholesky decomposition passes, if not it returns the rank at which the decomposition failed.

Perform Cholesky Decomposition

ArrayFire Docs

This function decomposes a positive definite matrix A into two triangular matrices.

choleskyInplace Source #

Arguments

:: AFType a 
=> Array a

is the input matrix on entry. It contains the triangular matrix on exit.

-> Bool

a boolean determining if in is upper or lower triangular

-> Int

is 0 if cholesky decomposition passes, if not it returns the rank at which the decomposition failed.

Perform Cholesky Decomposition

ArrayFire Docs

C Interface for in place cholesky decomposition.

solve Source #

Arguments

:: 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

ArrayFire Docs

solveLU Source #

Arguments

:: 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.

ArrayFire Docs

inverse Source #

Arguments

:: AFType a 
=> Array a

is input matrix

-> MatProp

determining various properties of matrix in

-> Array a

will contain the inverse of matrix in

Invert a matrix.

ArrayFire Docs

C Interface for inverting a matrix.

rank Source #

Arguments

:: 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

ArrayFire Docs

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

ArrayFire Docs

This function uses af::qr to find the rank of the input matrix within the given tolerance.

det Source #

Arguments

:: AFType a 
=> Array a

is input matrix

-> (Double, Double)

will contain the real and imaginary part of the determinant of in

Find the determinant of a Matrix

ArrayFire Docs

C Interface for finding the determinant of a matrix.

norm Source #

Arguments

:: AFType a 
=> Array a

is the input matrix

-> NormType

specifies the NormType

-> 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.

ArrayFire Docs

This function can return the norm using various metrics based on the type paramter.

isLAPACKAvailable Source #

Arguments

:: Bool

Returns if LAPACK is available

Is LAPACK available

ArrayFire Docs