Do more. Code less. Free software for GPU computing.
<scroll to top>
Functions

Factorizations: LU, QR, Cholesky, singular values, eigenvalues, Hessenberg

Linear Algebra

Functions

array lu (const array &in)
 LU factorization (packed)
void lu (array &lower, array &upper, const array &in)
 LU factorization.
void lu (array &lower, array &upper, array &pivot, const array &in)
 LU factorization (with pivoting)
array qr (const array &in)
 QR factorization (packed).
void qr (array &q, array &r, const array &in)
 QR factorization.
void qr (array &q, array &r, array &tau, const array &in)
 QR factorization with tau.
array cholesky (unsigned &info, const array &X, bool is_upper=true)
 Cholesky decomposition ("Y^T * Y == X").
array hessenberg (const array &in)
 Hessenberg matrix form.
void hessenberg (array &h, array &q, const array &in)
 Hessenberg matrix h with unitary permutation matrix q.
array eigen (const array &in, bool is_diag=false)
 Eigenvalues.
void eigen (array &values, array &vectors, const array &in)
 Eigenvalues and eigenvectors.
array svd (const array &in, bool is_diag=false)
 Singular values.
void svd (array &s, array &u, array &v, const array &in)
 Singular values with unitary bases: in = u * s * v.

Function Documentation

array af::lu ( const array &  in)

LU factorization (packed)

Parameters:
[in]in
Returns:
lower and upper matrix combined (getrf)
Examples:
examples/misc/lin_algebra.cpp.
void af::lu ( array &  lower,
array &  upper,
const array &  in 
)

LU factorization.

Parameters:
[out]lowertriangular matrix
[out]uppertriangular matrix
[in]in
void af::lu ( array &  lower,
array &  upper,
array &  pivot,
const array &  in 
)

LU factorization (with pivoting)

Parameters:
[out]lowertriangular matrix
[out]uppertriangular matrix
[out]pivotindices
[in]in
array af::qr ( const array &  in)

QR factorization (packed).

Parameters:
[in]in
Returns:
Q,R factors (geqrf)
void af::qr ( array &  q,
array &  r,
const array &  in 
)

QR factorization.

Parameters:
[out]qunitary rotation matrix
[out]rupper triangular
[in]in
void af::qr ( array &  q,
array &  r,
array &  tau,
const array &  in 
)

QR factorization with tau.

Parameters:
[out]qorthogonal matrix
[out]rupper triangular matrix
[out]tauAdditional information about q
[in]in
array af::cholesky ( unsigned &  info,
const array &  X,
bool  is_upper = true 
)

Cholesky decomposition ("Y^T * Y == X").

Parameters:
[out]infodetails on result of decomposition
[in]X
[in]is_uppertrue for upper triangular solution, false for lower triangular solution
Returns:
upper or lower triangular factor
array af::hessenberg ( const array &  in)

Hessenberg matrix form.

Parameters:
[in]in
Returns:
Hessenberg form
void af::hessenberg ( array &  h,
array &  q,
const array &  in 
)

Hessenberg matrix h with unitary permutation matrix q.

such that in = q * h * q^T. Requires ArrayFire Pro.

Parameters:
[out]hHessenberg matrix
[out]qunitary matrix
[in]in
array af::eigen ( const array &  in,
bool  is_diag = false 
)

Eigenvalues.

Requires ArrayFire Pro.

Parameters:
[in]in
[in]is_diagtrue then produce diagonal matrix of eigenvalues, false then vector of eigenvalues
Returns:
eigenvalues
Examples:
examples/misc/lin_algebra.cpp.
void af::eigen ( array &  values,
array &  vectors,
const array &  in 
)

Eigenvalues and eigenvectors.

Requires ArrayFire Pro.

Parameters:
[out]values
[out]vectors
[in]in
array af::svd ( const array &  in,
bool  is_diag = false 
)

Singular values.

Double-precision or complex values require ArrayFire Pro.

Parameters:
[in]in
[in]is_diagtrue then produce diagonal matrix of eigenvalues, false then vector of eigenvalues
Returns:
singular values
void af::svd ( array &  s,
array &  u,
array &  v,
const array &  in 
)

Singular values with unitary bases: in = u * s * v.

Parameters:
[out]ssingular values
[out]uleft unitary matrix
[out]vright unitary matrix
[in]inThe input matrix.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines