hasktorch-indef-0.0.1.0: Core Hasktorch abstractions wrapping FFI bindings

Copyright(c) Sam Stites 2017
LicenseBSD3
Maintainersam@stites.io
Stabilityexperimental
Portabilitynon-portable
Safe HaskellNone
LanguageHaskell2010

Torch.Indef.Static.Tensor.Math.Lapack

Description

 
Synopsis

Documentation

getri :: Tensor d -> Tensor d' Source #

Docs taken from MAGMA documentation at: http://icl.cs.utk.edu/projectsfiles/magma/doxygen/group__magma__getri.html

getri computes the inverse of a matrix using the LU factorization computed by getrf.

This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A).

Note that it is generally both faster and more accurate to use gesv, or getrf and getrs, to solve the system AX = B, rather than inverting the matrix and multiplying to form X = inv(A)*B. Only in special instances should an explicit inverse be computed with this routine.

getri_ :: Tensor d -> IO () Source #

inplace version of getri

potrf Source #

Arguments

:: Tensor d

matrix to decompose

-> Triangle

which triangle should be used.

-> Tensor d' 

Cholesky Decomposition of 2D tensor A. The matrix, A, has to be a positive-definite and either symmetric or complex Hermitian.

The factorization has the form

  A = U**T * U,   if UPLO = U, or
  A = L  * L**T,  if UPLO = L,

where U is an upper triangular matrix and L is lower triangular.

potrf_ Source #

Arguments

:: Tensor d

matrix to decompose

-> Triangle

which triangle should be used.

-> IO () 

infix version of potrf.

potri :: Tensor d -> Triangle -> Tensor d' Source #

Returns the inverse of 2D tensor A given its Cholesky decomposition chol.

Square matrix chol should be triangular.

Triangle specifies matrix chol as either upper or lower triangular.

potri_ :: Tensor d -> Triangle -> IO () Source #

inplace version of potri.

potrs Source #

Arguments

:: Tensor d

Tensor B

-> Tensor d'

Cholesky decomposition chol

-> Triangle

which triangle to use (upper or lower)

-> Tensor d'' 

Returns the solution to linear system AX = B using the Cholesky decomposition chol of 2D tensor A.

Square matrix chol should be triangular; and, righthand side matrix B should be of full rank.

Triangle specifies matrix chol as either upper or lower triangular.

potrs_ Source #

Arguments

:: Tensor d

Tensor B

-> Tensor d'

Cholesky decomposition chol

-> Triangle

which triangle to use (upper or lower)

-> IO () 

Inplace version of potri. Mutating tensor B in place.

qr :: Tensor d -> (Tensor d', Tensor d'') Source #

Compute a QR decomposition of the matrix x: matrices q and r such that x = q * r, with q orthogonal and r upper triangular. This returns the thin (reduced) QR factorization.

Note that precision may be lost if the magnitudes of the elements of x are large.

Note also that, while it should always give you a valid decomposition, it may not give you the same one across platforms - it will depend on your LAPACK implementation.

Note: Irrespective of the original strides, the returned matrix q will be transposed, i.e. with strides 1, m instead of m, 1.

qr_ :: (Tensor d, Tensor d') -> Tensor d'' -> IO () Source #

Inplace version of qr

_geqrf :: Tensor d -> Tensor d' -> Tensor d'' -> IO () Source #

This is a low-level function for calling LAPACK directly. You'll generally want to use qr instead.

Computes a QR decomposition of a, but without constructing Q and R as explicit separate matrices. Rather, this directly calls the underlying LAPACK function ?geqrf which produces a sequence of "elementary reflectors". See LAPACK documentation from MKL for further details and <http://icl.cs.utk.edu/projectsfiles/magma/doxygen/group__magma__geqrf.html for MAGMA documentation>.

Note that, because this is low-level code, hasktorch just calls Torch directly.

geev Source #

Arguments

:: Tensor d

square matrix to get eigen{values/vectors} of.

-> EigenReturn

whether or not to return eigenvectors.

-> (Tensor d', Maybe (Tensor d''))

(e, V) standing for eigenvalues and eigenvectors

(e, V) <- geev A returns eigenvalues and eigenvectors of a general real square matrix A.

A and V are m × m matrices and e is an m-dimensional vector.

This function calculates all right eigenvalues (and vectors) of A such that A = V diag(e) V.

The EigenReturn argument defines computation of eigenvectors or eigenvalues only. It determines if only eigenvalues are computed or if both eigenvalues and eigenvectors are computed.

The eigenvalues returned follow LAPACK convention and are returned as complex (real/imaginary) pairs of numbers (2 * m dimensional tensor).

Also called the "eig" fuction in torch.

geev_ Source #

Arguments

:: (Tensor d, Tensor d')

(e, V) standing for eigenvalues and eigenvectors

-> Tensor d''

square matrix to get eigen{values/vectors} of.

-> EigenReturn

whether or not to return eigenvectors.

-> IO () 

In-place version of geev.

Note: Irrespective of the original strides, the returned matrix V will be transposed, i.e. with strides 1, m instead of m, 1.

eig :: Tensor d -> EigenReturn -> (Tensor d', Maybe (Tensor d'')) Source #

alias to geev to match Torch naming conventions.

eig_ :: (Tensor d, Tensor d') -> Tensor d'' -> EigenReturn -> IO () Source #

alias to geev_ to match Torch naming conventions.

syev Source #

Arguments

:: Tensor d

square matrix to get eigen{values/vectors} of.

-> EigenReturn

whether or not to return eigenvectors.

-> Triangle

whether the upper or lower triangle should be used

-> (Tensor d', Maybe (Tensor d''))

(e, V) standing for eigenvalues and eigenvectors

(e, V) <- syev A returns eigenvalues and eigenvectors of a symmetric real matrix A.

A and V are m × m matrices and e is a m-dimensional vector.

This function calculates all eigenvalues (and vectors) of A such that A = V diag(e) V.

The EigenReturn argument defines computation of eigenvectors or eigenvalues only.

Since the input matrix A is supposed to be symmetric, only one triangular portion is used. The Triangle argument indicates if this should be the upper or lower triangle.

syev_ Source #

Arguments

:: (Tensor d, Tensor d')

(e, V) standing for eigenvalues and eigenvectors

-> Tensor d''

square matrix to get eigen{values/vectors} of.

-> EigenReturn

whether or not to return eigenvectors.

-> Triangle

whether the upper or lower triangle should be used

-> IO () 

Inplace version of syev

Note: Irrespective of the original strides, the returned matrix V will be transposed, i.e. with strides 1, m instead of m, 1.

symeig :: Tensor d -> EigenReturn -> Triangle -> (Tensor d', Maybe (Tensor d'')) Source #

alias to syev to match Torch naming conventions.

symeig_ :: (Tensor d, Tensor d') -> Tensor d'' -> EigenReturn -> Triangle -> IO () Source #

alias to syev to match Torch naming conventions.

gesv Source #

Arguments

:: Tensor '[m, k]
B
-> Tensor '[m, m]
A
-> (Tensor '[m, k], Tensor '[m, m])
(X, LU)

(X, LU) <- gesv B A returns the solution of AX = B and LU contains L and U factors for LU factorization of A.

A has to be a square and non-singular matrix (a 2D tensor). A and LU are m × m, X is m × k and B is m × k.

gesv_ Source #

Arguments

:: (Tensor '[m, k], Tensor '[m, m])
(X, LU)
-> Tensor '[m, k]
B
-> Tensor '[m, m]
A
-> IO () 

Inplace version of gesv.

In this case x and lu will be used for temporary storage and returning the result.

  • x will contain the solution X.
  • lu will contain L and U factors for LU factorization of A.

Note: Irrespective of the original strides, the returned matrices x and lu will be transposed, i.e. with strides 1, m instead of m, 1.

gels :: Tensor d -> Tensor d' -> (Tensor d'', Tensor d''') Source #

Solution of least squares and least norm problems for a full rank m × n matrix A.

  • If n ≤ m, then solve ||AX-B||_F.
  • If n > m, then solve min ||X||_F such that AX = B.

On return, first n rows of x matrix contains the solution and the rest contains residual information. Square root of sum squares of elements of each column of x starting at row n + 1 is the residual for corresponding column.

gels_ Source #

Arguments

:: (Tensor d, Tensor d')
(resb, resa)
-> Tensor d''
matrix b
-> Tensor d'''
matrix a
-> IO () 

Inplace version of gels.

Note: Irrespective of the original strides, the returned matrices resb and resa will be transposed, i.e. with strides 1, m instead of m, 1.

gesvd :: Tensor d -> ComputeSingularValues -> (Tensor d', Tensor d'', Tensor d''') Source #

(U, S, V) <- svd A returns the singular value decomposition of a real matrix A of size n × m such that A = USV'*.

U is n × n, S is n × m and V is m × m.

The ComputeSingularValues argument represents the number of singular values to be computed. SomeSVs stands for "some" (FIXME: figure out what that means) and AllSVs stands for all.

gesvd_ Source #

Arguments

:: (Tensor d, Tensor d', Tensor d'')

(u, s, v)

-> Tensor d'''

m

-> ComputeSingularValues

Whether to compute all or some of the singular values

-> IO () 

Inplace version of gesvd.

Note: Irrespective of the original strides, the returned matrix U will be transposed, i.e. with strides 1, n instead of n, 1.

gesvd2 Source #

Arguments

:: Tensor d

m

-> ComputeSingularValues

Whether to compute all or some of the singular values

-> (Tensor d', Tensor d'', Tensor d''', Tensor d'''')

(u, s, v, a)

gesvd, computing A = U*Σ*transpose(V).

NOTE: "gesvd, computing A = U*Σ*transpose(V)." is only inferred documentation. This documentation was made by stites, inferring from the description of the gesvd docs at <https://software.intel.com/en-us/mkl-developer-reference-c-gesvd the intel mkl documentation>.

gesvd2_ Source #

Arguments

:: (Tensor d, Tensor d', Tensor d'', Tensor d''')

(u, s, v, a)

-> Tensor d''''

m

-> ComputeSingularValues

Whether to compute all or some of the singular values

-> IO () 

Inplace version of gesvd2_.

data Triangle Source #

Argument to specify whether the upper or lower triangular decomposition should be used in potrf and potrf_.

Constructors

Upper

use upper triangular matrix

Lower

use lower triangular matrix

data EigenReturn Source #

Argument to be passed to geev, syev, and their inplace variants. Determines if the a function should only compute eigenvalues or both eigenvalues and eigenvectors.

data ComputeSingularValues Source #

Represents the number of singular values to be computed in gesvd and gesvd2. While fairly opaque about how many values are computed, Torch says we either compute "some" or all of the values.

Constructors

SomeSVs 
AllSVs