MAGMA  2.7.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
sy/heevdx: Solves using divide-and-conquer (expert)

Functions

magma_int_t magma_cheevdx (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *w, magmaFloatComplex *work, magma_int_t lwork, float *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 CHEEVDX computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_cheevdx_2stage (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *W, magmaFloatComplex *work, magma_int_t lwork, float *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 CHEEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_cheevdx_2stage_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *W, magmaFloatComplex *work, magma_int_t lwork, float *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 CHEEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_cheevdx_gpu (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex_ptr dA, magma_int_t ldda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, float *w, magmaFloatComplex *wA, magma_int_t ldwa, magmaFloatComplex *work, magma_int_t lwork, float *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 CHEEVDX_GPU computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_cheevdx_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *w, magmaFloatComplex *work, magma_int_t lwork, float *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 CHEEVD computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_dsyevdx (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *w, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 DSYEVDX computes selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_dsyevdx_2stage (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *W, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 DSYEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_dsyevdx_2stage_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *W, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 DSYEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_dsyevdx_gpu (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, double *w, double *wA, magma_int_t ldwa, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 DSYEVDX computes selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_dsyevdx_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *w, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 DSYEVD computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_ssyevdx (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *w, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 SSYEVDX computes selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_ssyevdx_2stage (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *W, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 SSYEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_ssyevdx_2stage_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *W, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 SSYEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_ssyevdx_gpu (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloat_ptr dA, magma_int_t ldda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, float *w, float *wA, magma_int_t ldwa, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 SSYEVDX computes selected eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_ssyevdx_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *m, float *w, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 SSYEVD computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A. More...
 
magma_int_t magma_zheevdx (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *w, magmaDoubleComplex *work, magma_int_t lwork, double *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 ZHEEVDX computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_zheevdx_2stage (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *W, magmaDoubleComplex *work, magma_int_t lwork, double *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 ZHEEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_zheevdx_2stage_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *W, magmaDoubleComplex *work, magma_int_t lwork, double *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 ZHEEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_zheevdx_gpu (magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex_ptr dA, magma_int_t ldda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, double *w, magmaDoubleComplex *wA, magma_int_t ldwa, magmaDoubleComplex *work, magma_int_t lwork, double *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 ZHEEVDX_GPU computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 
magma_int_t magma_zheevdx_m (magma_int_t ngpu, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *m, double *w, magmaDoubleComplex *work, magma_int_t lwork, double *rwork, magma_int_t lrwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info)
 ZHEEVD computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A. More...
 

Detailed Description

Function Documentation

magma_int_t magma_cheevdx ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  w,
magmaFloatComplex *  work,
magma_int_t  lwork,
float *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

CHEEVDX computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wREAL array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) COMPLEX array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( N + N*NB, 2*N + N**2 ). NB can be obtained through magma_get_chetrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) REAL array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_cheevdx_2stage ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  W,
magmaFloatComplex *  work,
magma_int_t  lwork,
float *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

CHEEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WREAL array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) COMPLEX array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_cbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) REAL array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_cheevdx_2stage_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  W,
magmaFloatComplex *  work,
magma_int_t  lwork,
float *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

CHEEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WREAL array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) COMPLEX array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_cbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) REAL array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_cheevdx_gpu ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex_ptr  dA,
magma_int_t  ldda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  mout,
float *  w,
magmaFloatComplex *  wA,
magma_int_t  ldwa,
magmaFloatComplex *  work,
magma_int_t  lwork,
float *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

CHEEVDX_GPU computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dACOMPLEX array on the GPU, dimension (LDDA, N). On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first mout columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]lddaINTEGER The leading dimension of the array DA. LDDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]moutINTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1.
[out]wREAL array, dimension (N) If INFO = 0, the required mout eigenvalues in ascending order.
wA(workspace) COMPLEX array, dimension (LDWA, N)
[in]ldwaINTEGER The leading dimension of the array wA. LDWA >= max(1,N).
[out]work(workspace) COMPLEX array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( N + N*NB, 2*N + N**2 ). NB can be obtained through magma_get_chetrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) REAL array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_cheevdx_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaFloatComplex *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  w,
magmaFloatComplex *  work,
magma_int_t  lwork,
float *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

CHEEVD computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, A contains the orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wREAL array, dimension (N) If INFO = 0, the eigenvalues in ascending order.
[out]work(workspace) COMPLEX array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( N + N*NB, 2*N + N**2 ). NB can be obtained through magma_get_chetrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) REAL array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_dsyevdx ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
double *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  w,
double *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

DSYEVDX 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. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ADOUBLE PRECISION array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, A contains the orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. If N > 0: 1 <= IL <= IU <= N. Otherwise: IL = 1 and IU = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wDOUBLE PRECISION array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= 2*N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( 2*N + N*NB, 1 + 6*N + 2*N**2 ). NB can be obtained through magma_get_dsytrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_dsyevdx_2stage ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
double *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  W,
double *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

DSYEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ADOUBLE PRECISION array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WDOUBLE PRECISION array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_dbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_dsyevdx_2stage_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
double *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  W,
double *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

DSYEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ADOUBLE PRECISION array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WDOUBLE PRECISION array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_dbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_dsyevdx_gpu ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaDouble_ptr  dA,
magma_int_t  ldda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  mout,
double *  w,
double *  wA,
magma_int_t  ldwa,
double *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

DSYEVDX 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. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dADOUBLE PRECISION array on the GPU, dimension (LDDA, N). On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first mout columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]lddaINTEGER The leading dimension of the array DA. LDDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]moutINTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1.
[out]wDOUBLE PRECISION array, dimension (N) If INFO = 0, the required mout eigenvalues in ascending order.
wA(workspace) DOUBLE PRECISION array, dimension (LDWA, N)
[in]ldwaINTEGER The leading dimension of the array wA. LDWA >= max(1,N).
[out]work(workspace) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= 2*N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( 2*N + N*NB, 1 + 6*N + 2*N**2 ). NB can be obtained through magma_get_dsytrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_dsyevdx_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
double *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  w,
double *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

DSYEVD computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ADOUBLE PRECISION array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, A contains the orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wDOUBLE PRECISION array, dimension (N) If INFO = 0, the eigenvalues in ascending order.
[out]work(workspace) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= 2*N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( 2*N + N*NB, 1 + 6*N + 2*N**2 ). NB can be obtained through magma_get_dsytrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_ssyevdx ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
float *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  w,
float *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

SSYEVDX 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. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]AREAL array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, A contains the orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. If N > 0: 1 <= IL <= IU <= N. Otherwise: IL = 1 and IU = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wREAL array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) REAL array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= 2*N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( 2*N + N*NB, 1 + 6*N + 2*N**2 ). NB can be obtained through magma_get_ssytrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_ssyevdx_2stage ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
float *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  W,
float *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

SSYEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]AREAL array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WREAL array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) REAL array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_sbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_ssyevdx_2stage_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
float *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  W,
float *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

SSYEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]AREAL array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WREAL array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) REAL array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_sbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_ssyevdx_gpu ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaFloat_ptr  dA,
magma_int_t  ldda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  mout,
float *  w,
float *  wA,
magma_int_t  ldwa,
float *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

SSYEVDX 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. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dAREAL array on the GPU, dimension (LDDA, N). On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first mout columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]lddaINTEGER The leading dimension of the array DA. LDDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]moutINTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1.
[out]wREAL array, dimension (N) If INFO = 0, the required mout eigenvalues in ascending order.
wA(workspace) REAL array, dimension (LDWA, N)
[in]ldwaINTEGER The leading dimension of the array wA. LDWA >= max(1,N).
[out]work(workspace) REAL array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= 2*N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( 2*N + N*NB, 1 + 6*N + 2*N**2 ). NB can be obtained through magma_get_ssytrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK and IWORK arrays, returns these values as the first entries of the WORK and IWORK arrays, and no error message related to LWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_ssyevdx_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
float *  A,
magma_int_t  lda,
float  vl,
float  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
float *  w,
float *  work,
magma_int_t  lwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

SSYEVD computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]AREAL array, dimension (LDA, N) On entry, the symmetric matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, A contains the orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlREAL
[in]vuREAL If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wREAL array, dimension (N) If INFO = 0, the eigenvalues in ascending order.
[out]work(workspace) REAL array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= 2*N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( 2*N + N*NB, 1 + 6*N + 2*N**2 ). NB can be obtained through magma_get_ssytrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_zheevdx ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  w,
magmaDoubleComplex *  work,
magma_int_t  lwork,
double *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

ZHEEVDX computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX_16 array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wDOUBLE PRECISION array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) COMPLEX_16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( N + N*NB, 2*N + N**2 ). NB can be obtained through magma_get_zhetrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) DOUBLE PRECISION array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_zheevdx_2stage ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  W,
magmaDoubleComplex *  work,
magma_int_t  lwork,
double *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

ZHEEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX_16 array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WDOUBLE PRECISION array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) COMPLEX_16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_zbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) DOUBLE PRECISION array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_zheevdx_2stage_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  W,
magmaDoubleComplex *  work,
magma_int_t  lwork,
double *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

ZHEEVD_2STAGE_M computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

It uses a two-stage algorithm for the tridiagonalization. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX_16 array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first m columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]WDOUBLE PRECISION array, dimension (N) If INFO = 0, the required m eigenvalues in ascending order.
[out]work(workspace) COMPLEX_16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= LWSTG2 + N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= LWSTG2 + 2*N + N**2. where LWSTG2 is the size needed to store the matrices of stage 2 and is returned by magma_zbulge_getlwstg2.
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) DOUBLE PRECISION array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_zheevdx_gpu ( magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex_ptr  dA,
magma_int_t  ldda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  mout,
double *  w,
magmaDoubleComplex *  wA,
magma_int_t  ldwa,
magmaDoubleComplex *  work,
magma_int_t  lwork,
double *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

ZHEEVDX_GPU computes selected eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues. If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]dACOMPLEX_16 array on the GPU, dimension (LDDA, N). On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, the first mout columns of A contains the required orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]lddaINTEGER The leading dimension of the array DA. LDDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]moutINTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1.
[out]wDOUBLE PRECISION array, dimension (N) If INFO = 0, the required mout eigenvalues in ascending order.
wA(workspace) COMPLEX_16 array, dimension (LDWA, N)
[in]ldwaINTEGER The leading dimension of the array wA. LDWA >= max(1,N).
[out]work(workspace) COMPLEX_16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( N + N*NB, 2*N + N**2 ). NB can be obtained through magma_get_zhetrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) DOUBLE PRECISION array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.

magma_int_t magma_zheevdx_m ( magma_int_t  ngpu,
magma_vec_t  jobz,
magma_range_t  range,
magma_uplo_t  uplo,
magma_int_t  n,
magmaDoubleComplex *  A,
magma_int_t  lda,
double  vl,
double  vu,
magma_int_t  il,
magma_int_t  iu,
magma_int_t *  m,
double *  w,
magmaDoubleComplex *  work,
magma_int_t  lwork,
double *  rwork,
magma_int_t  lrwork,
magma_int_t *  iwork,
magma_int_t  liwork,
magma_int_t *  info 
)

ZHEEVD computes all eigenvalues and, optionally, eigenvectors of a complex Hermitian matrix A.

If eigenvectors are desired, it uses a divide and conquer algorithm.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Parameters
[in]ngpuINTEGER Number of GPUs to use. ngpu > 0.
[in]jobzmagma_vec_t
  • = MagmaNoVec: Compute eigenvalues only;
  • = MagmaVec: Compute eigenvalues and eigenvectors.
[in]rangemagma_range_t
  • = MagmaRangeAll: all eigenvalues will be found.
  • = MagmaRangeV: all eigenvalues in the half-open interval (VL,VU] will be found.
  • = MagmaRangeI: the IL-th through IU-th eigenvalues will be found.
[in]uplomagma_uplo_t
  • = MagmaUpper: Upper triangle of A is stored;
  • = MagmaLower: Lower triangle of A is stored.
[in]nINTEGER The order of the matrix A. N >= 0.
[in,out]ACOMPLEX_16 array, dimension (LDA, N) On entry, the Hermitian matrix A. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A. If UPLO = MagmaLower, the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A. On exit, if JOBZ = MagmaVec, then if INFO = 0, A contains the orthonormal eigenvectors of the matrix A. If JOBZ = MagmaNoVec, then on exit the lower triangle (if UPLO=MagmaLower) or the upper triangle (if UPLO=MagmaUpper) of A, including the diagonal, is destroyed.
[in]ldaINTEGER The leading dimension of the array A. LDA >= max(1,N).
[in]vlDOUBLE PRECISION
[in]vuDOUBLE PRECISION If RANGE=MagmaRangeV, the lower and upper bounds of the interval to be searched for eigenvalues. VL < VU. Not referenced if RANGE = MagmaRangeAll or MagmaRangeI.
[in]ilINTEGER
[in]iuINTEGER If RANGE=MagmaRangeI, the indices (in ascending order) of the smallest and largest eigenvalues to be returned. 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if RANGE = MagmaRangeAll or MagmaRangeV.
[out]mINTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1.
[out]wDOUBLE PRECISION array, dimension (N) If INFO = 0, the eigenvalues in ascending order.
[out]work(workspace) COMPLEX_16 array, dimension (MAX(1,LWORK)) On exit, if INFO = 0, WORK[0] returns the optimal LWORK.
[in]lworkINTEGER The length of the array WORK.
  • If N <= 1, LWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LWORK >= N + N*NB.
  • If JOBZ = MagmaVec and N > 1, LWORK >= max( N + N*NB, 2*N + N**2 ). NB can be obtained through magma_get_zhetrd_nb(N).
    If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]rwork(workspace) DOUBLE PRECISION array, dimension (LRWORK) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK.
[in]lrworkINTEGER The dimension of the array RWORK.
  • If N <= 1, LRWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LRWORK >= N.
  • If JOBZ = MagmaVec and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
    If LRWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]iwork(workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK.
[in]liworkINTEGER The dimension of the array IWORK.
  • If N <= 1, LIWORK >= 1.
  • If JOBZ = MagmaNoVec and N > 1, LIWORK >= 1.
  • If JOBZ = MagmaVec and N > 1, LIWORK >= 3 + 5*N.
    If LIWORK = -1, then a workspace query is assumed; the routine only calculates the optimal sizes of the WORK, RWORK and IWORK arrays, returns these values as the first entries of the WORK, RWORK and IWORK arrays, and no error message related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]infoINTEGER
  • = 0: successful exit
  • < 0: if INFO = -i, the i-th argument had an illegal value
  • > 0: if INFO = i and JOBZ = MagmaNoVec, then the algorithm failed to converge; i off-diagonal elements of an intermediate tridiagonal form did not converge to zero; if INFO = i and JOBZ = MagmaVec, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns INFO/(N+1) through mod(INFO,N+1).

Further Details

Based on contributions by Jeff Rutter, Computer Science Division, University of California at Berkeley, USA

Modified description of INFO. Sven, 16 Feb 05.