mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Templatized Lapack interface

Templatized interface to the Lapack library

Functions

template<typename dataT >
dataT mx::math::lamch (char CMACH)
 Determine machine parameters. More...
 
template<typename dataT >
MXLAPACK_INT mx::math::sytrd (char UPLO, MXLAPACK_INT N, dataT *A, MXLAPACK_INT LDA, dataT *D, dataT *E, dataT *TAU, dataT *WORK, MXLAPACK_INT LWORK, MXLAPACK_INT INFO)
 Reduce a real symmetric matrix to real symmetric tridiagonal form by an orthogonal similarity transformation. More...
 
template<typename dataT >
MXLAPACK_INT mx::math::syevr (char JOBZ, char RANGE, char UPLO, MXLAPACK_INT N, dataT *A, MXLAPACK_INT LDA, dataT VL, dataT VU, MXLAPACK_INT IL, MXLAPACK_INT IU, dataT ABSTOL, MXLAPACK_INT *M, dataT *W, dataT *Z, MXLAPACK_INT LDZ, MXLAPACK_INT *ISUPPZ, dataT *WORK, MXLAPACK_INT LWORK, MXLAPACK_INT *IWORK, MXLAPACK_INT LIWORK)
 Compute selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix. More...
 
template<typename dataT >
MXLAPACK_INT mx::math::gesvd (char JOBU, char JOBVT, MXLAPACK_INT M, MXLAPACK_INT N, dataT *A, MXLAPACK_INT LDA, dataT *S, dataT *U, MXLAPACK_INT LDU, dataT *VT, MXLAPACK_INT LDVT, dataT *WORK, MXLAPACK_INT LWORK)
 Compute the singular value decomposition (SVD) of a real matrix. More...
 

Function Documentation

◆ gesvd()

template<typename dataT >
MXLAPACK_INT mx::math::gesvd ( char  JOBU,
char  JOBVT,
MXLAPACK_INT  M,
MXLAPACK_INT  N,
dataT *  A,
MXLAPACK_INT  LDA,
dataT *  S,
dataT *  U,
MXLAPACK_INT  LDU,
dataT *  VT,
MXLAPACK_INT  LDVT,
dataT *  WORK,
MXLAPACK_INT  LWORK 
)

Compute the singular value decomposition (SVD) of a real matrix.

xGESVD computes the singular value decomposition (SVD) of a real M-by-N matrix A, optionally computing the left and/or right singular vectors. The SVD is written

\( A = U * \Sigma * V^T \)

where \( \Sigma \) is an M-by-N matrix which is zero except for its min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and V is an N-by-N orthogonal matrix. The diagonal elements of \( \Sigma \) are the singular values of A; they are real and non-negative, and are returned in descending order. The first min(m,n) columns of U and V are the left and right singular vectors of A.

Note that the routine returns \( V^T \), not \( V \).

See more details: http://www.netlib.org/lapack/explore-html/d8/d49/sgesvd_8f.html. This documentation is taken from there.

Template Parameters
dataTis the data type of the arrays, and determines which underlying Lapack routine is called.
Parameters
[in]JOBU

(char) Specifies options for computing all or part of the matrix U:
= 'A': all M columns of U are returned in array U
= 'S': the first min(m,n) columns of U (the left singular vectors) are returned in the array U
= 'O': the first min(m,n) columns of U (the left singular
vectors) are overwritten on the array A
= 'N': no columns of U (no left singular vectors) are computed.

[in]JOBVT

(char) Specifies options for computing all or part of the matrix V**T: = 'A': all N rows of V**T are returned in the array VT
= 'S': the first min(m,n) rows of V**T (the right singular vectors) are returned in the array VT
= 'O': the first min(m,n) rows of V**T (the right singular vectors) are overwritten on the array A
= 'N': no rows of V**T (no right singular vectors) are computed.

JOBVT and JOBU cannot both be 'O'.

[in]M

(MXLAPACK_INT) The number of rows of the input matrix A. M >= 0.

[in]N

(MXLAPACK_INT) The number of columns of the input matrix A. N >= 0.

[in,out]A

(dataT *, dimension (LDA,N)) On entry: the M-by-N matrix A.
On exit:
if JOBU = 'O', A is overwritten with the first min(m,n) columns of U (the left singular vectors, stored columnwise)
if JOBVT = 'O', A is overwritten with the first min(m,n) rows of V**T (the right singular vectors, stored rowwise)
if JOBU != 'O' and JOBVT .ne. 'O', the contents of A are destroyed.

[in]LDA

(MXLAPACK_INT) The leading dimension of the array A. LDA >= max(1,M).

[out]S

(dataT *, dimension (min(M,N)) The singular values of A, sorted so that S(i) >= S(i+1).

[out]U

(dataT *, dimension (LDU,UCOL)) (LDU,M) if JOBU = 'A' or (LDU,min(M,N)) if JOBU = 'S'.
If JOBU = 'A', U contains the M-by-M orthogonal matrix U
if JOBU = 'S', U contains the first min(m,n) columns of U (the left singular vectors, stored columnwise)
if JOBU = 'N' or 'O', U is not referenced.

[in]LDU

(MXLAPACK_INT) The leading dimension of the array U.
LDU >= 1
if JOBU == 'S' or 'A', LDU >= M.

[out]VT

(dataT *, dimension (LDVT,N)) If JOBVT = 'A', VT contains the N-by-N orthogonal matrix V**T
if JOBVT = 'S', VT contains the first min(m,n) rows of V**T (the right singular vectors, stored rowwise)
if JOBVT = 'N' or 'O', VT is not referenced.

[in]LDVT

(MXLAPACK_INT) The leading dimension of the array VT.
LDVT >= 1
if JOBVT = 'A', LDVT >= N
if JOBVT = 'S', LDVT >= min(M,N).

[out]WORK

(dataT *, dimension (MAX(1,LWORK)) )
On exit, if INFO = 0, WORK[1] returns the optimal LWORK
if INFO > 0, WORK(2:MIN(M,N)) contains the unconverged superdiagonal elements of an upper bidiagonal matrix B whose diagonal is in S (not necessarily sorted). B satisfies A = U * B * VT, so it has the same singular values as A, and singular vectors related by U and VT.

[in]LWORK

(MXLAPACK_INT) The dimension of the array WORK.
LWORK >= MAX(1,5*MIN(M,N)) for the paths (see comments inside code):

  • PATH 1 (M much larger than N, JOBU='N')
  • PATH 1t (N much larger than M, JOBVT='N') LWORK >= MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N)) for the other paths For good performance, LWORK should generally be larger.

If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.

Returns
=0: successful exit.
<0: if INFO = -i, the i-th argument had an illegal value.
>0: if SBDSQR did not converge, INFO specifies how many superdiagonals of an MXLAPACK_INTermediate bidiagonal form B did not converge to zero. See the description of WORK above for details.

Definition at line 412 of file templateLapack.hpp.

◆ lamch()

template<typename dataT >
dataT mx::math::lamch ( char  CMACH)

Determine machine parameters.

Wrapper for Lapack xLAMCH

See more details at http://www.netlib.org/lapack/lapack-3.1.1/html/slamch.f.html.

  • Template Parameters
    dataTis the data type of the arrays, and determines which underlying Lapack routine is called.
    Parameters
    [in]CMACHSpecifies the value to be returned by SLAMCH:
        = 'E' or 'e',   returns eps, the relative machine precision  
        = 'S' or 's ,   returns sfmin, the safe minimum, such that 1/sfmin does not overflow  
        = 'B' or 'b',   returns base, the base of the machine 
        = 'P' or 'p',   returns eps*base  
        = 'N' or 'n',   returns t, the number of (base) digits in the mantissa  
        = 'R' or 'r',   returns rnd, 1.0 when rounding occurs in addition, 0.0 otherwise  
        = 'M' or 'm',   returns emin, the minimum exponent before (gradual) underflow  
        = 'U' or 'u',   returns rmin, the underflow threshold - base^(emin-1)  
        = 'L' or 'l',   returns emax, the largest exponent before overflow  
        = 'O' or 'o',   returns rmax, the overflow threshold  - (base^emax)*(1-eps)  
    
    Returns
    the value of the specified machine parameters for the specified precision

Definition at line 110 of file templateLapack.hpp.

◆ syevr()

template<typename dataT >
MXLAPACK_INT mx::math::syevr ( char  JOBZ,
char  RANGE,
char  UPLO,
MXLAPACK_INT  N,
dataT *  A,
MXLAPACK_INT  LDA,
dataT  VL,
dataT  VU,
MXLAPACK_INT  IL,
MXLAPACK_INT  IU,
dataT  ABSTOL,
MXLAPACK_INT *  M,
dataT *  W,
dataT *  Z,
MXLAPACK_INT  LDZ,
MXLAPACK_INT *  ISUPPZ,
dataT *  WORK,
MXLAPACK_INT  LWORK,
MXLAPACK_INT *  IWORK,
MXLAPACK_INT  LIWORK 
)

Compute selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix.

xSYEVR computes selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.

See more details: http://www.netlib.org/lapack/lapack-3.1.1/html/ssyevr.f.html.

Definition at line 223 of file templateLapack.hpp.

◆ sytrd()

template<typename dataT >
MXLAPACK_INT mx::math::sytrd ( char  UPLO,
MXLAPACK_INT  N,
dataT *  A,
MXLAPACK_INT  LDA,
dataT *  D,
dataT *  E,
dataT *  TAU,
dataT *  WORK,
MXLAPACK_INT  LWORK,
MXLAPACK_INT  INFO 
)

Reduce a real symmetric matrix to real symmetric tridiagonal form by an orthogonal similarity transformation.

xSYTRD reduces a real symmetric matrix A to real symmetric tridiagonal form T by an orthogonal similarity transformation:

\( Q^T * A * Q = T. \)

For more see: http://www.netlib.org/lapack/lapack-3.1.1/html/ssytrd.f.html

Template Parameters
dataTis the data type
Parameters
UPLO'U': Upper triangle of A is stored, 'L': Lower triangle of A is stored.
NThe order of the matrix A. N >= 0.
Aarray, dimension (LDA,N)
LDAThe leading dimension of the array A. LDA >= max(1,N).
D(output) array, dimension (N), the diagonal elements of the tridiagonal matrix T:
E(output) array, dimension (N-1), the off-diagonal elements of the tridiagonal matrix T:
TAU(output) array, dimension (N-1), the scalar factors of the elementary reflectors (see Further Details).
WORK(workspace/output) array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
LWORK(input) The dimension of the array WORK. LWORK >= 1. For optimum performance LWORK >= N*NB, where NB is the optimal blocksize.
INFO(output) 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value
Returns
the value of INFO from the LAPACK routine

Definition at line 183 of file templateLapack.hpp.