![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
magma_int_t | magma_chegvdx (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, magmaFloatComplex *B, magma_int_t ldb, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, 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) |
CHEGVDX computes selected eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_chegvdx_2stage (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, magmaFloatComplex *B, magma_int_t ldb, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, 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) |
CHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_chegvdx_2stage_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, magmaFloatComplex *B, magma_int_t ldb, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, 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) |
CHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_chegvdx_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaFloatComplex *A, magma_int_t lda, magmaFloatComplex *B, magma_int_t ldb, 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) |
CHEGVD computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_dsygvdx (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double *B, magma_int_t ldb, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, double *w, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info) |
DSYGVDX computes selected eigenvalues and, optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_dsygvdx_2stage (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double *B, magma_int_t ldb, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, double *w, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info) |
DSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_dsygvdx_2stage_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double *B, magma_int_t ldb, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, double *w, double *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info) |
DSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_dsygvdx_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, double *A, magma_int_t lda, double *B, magma_int_t ldb, 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) |
DSYGVD computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_ssygvdx (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float *B, magma_int_t ldb, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, float *w, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info) |
SSYGVDX computes selected eigenvalues and, optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_ssygvdx_2stage (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float *B, magma_int_t ldb, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, float *w, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info) |
SSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_ssygvdx_2stage_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float *B, magma_int_t ldb, float vl, float vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, float *w, float *work, magma_int_t lwork, magma_int_t *iwork, magma_int_t liwork, magma_int_t *info) |
SSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_ssygvdx_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, float *A, magma_int_t lda, float *B, magma_int_t ldb, 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) |
SSYGVD computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_zhegvdx (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, magmaDoubleComplex *B, magma_int_t ldb, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, 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) |
ZHEGVDX computes selected eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_zhegvdx_2stage (magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, magmaDoubleComplex *B, magma_int_t ldb, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, 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) |
ZHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_zhegvdx_2stage_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, magmaDoubleComplex *B, magma_int_t ldb, double vl, double vu, magma_int_t il, magma_int_t iu, magma_int_t *mout, 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) |
ZHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t | magma_zhegvdx_m (magma_int_t ngpu, magma_int_t itype, magma_vec_t jobz, magma_range_t range, magma_uplo_t uplo, magma_int_t n, magmaDoubleComplex *A, magma_int_t lda, magmaDoubleComplex *B, magma_int_t ldb, 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) |
ZHEGVD computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. | |
magma_int_t magma_chegvdx | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaFloatComplex * | A, | ||
magma_int_t | lda, | ||
magmaFloatComplex * | B, | ||
magma_int_t | ldb, | ||
float | vl, | ||
float | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
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 ) |
CHEGVDX computes selected eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved: = 1: A*x = (lambda)*B*x = 2: A*B*x = (lambda)*x = 3: B*A*x = (lambda)*x |
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**H*B*Z = I; if ITYPE = 3, Z**H*inv(B)*Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) REAL array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if CHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_chegvdx_2stage | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaFloatComplex * | A, | ||
magma_int_t | lda, | ||
magmaFloatComplex * | B, | ||
magma_int_t | ldb, | ||
float | vl, | ||
float | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
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 ) |
CHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) REAL array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if CHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_chegvdx_2stage_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaFloatComplex * | A, | ||
magma_int_t | lda, | ||
magmaFloatComplex * | B, | ||
magma_int_t | ldb, | ||
float | vl, | ||
float | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
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 ) |
CHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) REAL array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if CHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_chegvdx_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaFloatComplex * | A, | ||
magma_int_t | lda, | ||
magmaFloatComplex * | B, | ||
magma_int_t | ldb, | ||
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 ) |
CHEGVD computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**H*B*Z = I; if ITYPE = 3, Z**H*inv(B)*Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | m | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) REAL array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if CHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_dsygvdx | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
double * | A, | ||
magma_int_t | lda, | ||
double * | B, | ||
magma_int_t | ldb, | ||
double | vl, | ||
double | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
double * | w, | ||
double * | work, | ||
magma_int_t | lwork, | ||
magma_int_t * | iwork, | ||
magma_int_t | liwork, | ||
magma_int_t * | info ) |
DSYGVDX computes selected eigenvalues and, optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved: = 1: A*x = (lambda)*B*x = 2: A*B*x = (lambda)*x = 3: B*A*x = (lambda)*x |
[in] | range | magma_range_t
|
[in] | jobz | magma_vec_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | DOUBLE 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**T * B * Z = I; if ITYPE = 3, Z**T * inv(B) * Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | DOUBLE PRECISION array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**T * U or B = L * L**T. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if DSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_dsygvdx_2stage | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
double * | A, | ||
magma_int_t | lda, | ||
double * | B, | ||
magma_int_t | ldb, | ||
double | vl, | ||
double | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
double * | w, | ||
double * | work, | ||
magma_int_t | lwork, | ||
magma_int_t * | iwork, | ||
magma_int_t | liwork, | ||
magma_int_t * | info ) |
DSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | DOUBLE 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | DOUBLE PRECISION array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if DSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_dsygvdx_2stage_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
double * | A, | ||
magma_int_t | lda, | ||
double * | B, | ||
magma_int_t | ldb, | ||
double | vl, | ||
double | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
double * | w, | ||
double * | work, | ||
magma_int_t | lwork, | ||
magma_int_t * | iwork, | ||
magma_int_t | liwork, | ||
magma_int_t * | info ) |
DSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | DOUBLE 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | DOUBLE PRECISION array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if DSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_dsygvdx_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
double * | A, | ||
magma_int_t | lda, | ||
double * | B, | ||
magma_int_t | ldb, | ||
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 ) |
DSYGVD computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved: = 1: A*x = (lambda)*B*x = 2: A*B*x = (lambda)*x = 3: B*A*x = (lambda)*x |
[in] | range | magma_range_t
|
[in] | jobz | magma_vec_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | DOUBLE 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**T*B*Z = I; if ITYPE = 3, Z**T*inv(B)*Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | DOUBLE PRECISION array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**T*U or B = L*L**T. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | m | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if DSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_ssygvdx | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
float * | A, | ||
magma_int_t | lda, | ||
float * | B, | ||
magma_int_t | ldb, | ||
float | vl, | ||
float | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
float * | w, | ||
float * | work, | ||
magma_int_t | lwork, | ||
magma_int_t * | iwork, | ||
magma_int_t | liwork, | ||
magma_int_t * | info ) |
SSYGVDX computes selected eigenvalues and, optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved: = 1: A*x = (lambda)*B*x = 2: A*B*x = (lambda)*x = 3: B*A*x = (lambda)*x |
[in] | range | magma_range_t
|
[in] | jobz | magma_vec_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | REAL 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**T * B * Z = I; if ITYPE = 3, Z**T * inv(B) * Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | REAL array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**T * U or B = L * L**T. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if SSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_ssygvdx_2stage | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
float * | A, | ||
magma_int_t | lda, | ||
float * | B, | ||
magma_int_t | ldb, | ||
float | vl, | ||
float | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
float * | w, | ||
float * | work, | ||
magma_int_t | lwork, | ||
magma_int_t * | iwork, | ||
magma_int_t | liwork, | ||
magma_int_t * | info ) |
SSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | REAL 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | REAL array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if SSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_ssygvdx_2stage_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
float * | A, | ||
magma_int_t | lda, | ||
float * | B, | ||
magma_int_t | ldb, | ||
float | vl, | ||
float | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
float * | w, | ||
float * | work, | ||
magma_int_t | lwork, | ||
magma_int_t * | iwork, | ||
magma_int_t | liwork, | ||
magma_int_t * | info ) |
SSYGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | REAL 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | REAL array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if SSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_ssygvdx_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
float * | A, | ||
magma_int_t | lda, | ||
float * | B, | ||
magma_int_t | ldb, | ||
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 ) |
SSYGVD computes all the eigenvalues, and optionally, the eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be symmetric and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved: = 1: A*x = (lambda)*B*x = 2: A*B*x = (lambda)*x = 3: B*A*x = (lambda)*x |
[in] | range | magma_range_t
|
[in] | jobz | magma_vec_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | REAL 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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**T*B*Z = I; if ITYPE = 3, Z**T*inv(B)*Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | REAL array, dimension (LDB, N) On entry, the symmetric matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**T*U or B = L*L**T. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | REAL |
[in] | vu | REAL 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] | il | INTEGER |
[in] | iu | INTEGER 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] | m | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | REAL 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] | lwork | INTEGER The length of the array WORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if SSYEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_zhegvdx | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaDoubleComplex * | A, | ||
magma_int_t | lda, | ||
magmaDoubleComplex * | B, | ||
magma_int_t | ldb, | ||
double | vl, | ||
double | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
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 ) |
ZHEGVDX computes selected eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved: = 1: A*x = (lambda)*B*x = 2: A*B*x = (lambda)*x = 3: B*A*x = (lambda)*x |
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX_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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**H*B*Z = I; if ITYPE = 3, Z**H*inv(B)*Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX_16 array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= MOUT <= N. If RANGE = MagmaRangeAll, MOUT = N, and if RANGE = MagmaRangeI, MOUT = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) DOUBLE PRECISION array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if ZHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_zhegvdx_2stage | ( | magma_int_t | itype, |
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaDoubleComplex * | A, | ||
magma_int_t | lda, | ||
magmaDoubleComplex * | B, | ||
magma_int_t | ldb, | ||
double | vl, | ||
double | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
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 ) |
ZHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX_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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX_16 array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) DOUBLE PRECISION array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if ZHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_zhegvdx_2stage_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaDoubleComplex * | A, | ||
magma_int_t | lda, | ||
magmaDoubleComplex * | B, | ||
magma_int_t | ldb, | ||
double | vl, | ||
double | vu, | ||
magma_int_t | il, | ||
magma_int_t | iu, | ||
magma_int_t * | mout, | ||
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 ) |
ZHEGVDX_2STAGE computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX_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 matrix Z of eigenvectors. The eigenvectors are normalized as follows:
|
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX_16 array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | mout | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) DOUBLE PRECISION array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if ZHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.
magma_int_t magma_zhegvdx_m | ( | magma_int_t | ngpu, |
magma_int_t | itype, | ||
magma_vec_t | jobz, | ||
magma_range_t | range, | ||
magma_uplo_t | uplo, | ||
magma_int_t | n, | ||
magmaDoubleComplex * | A, | ||
magma_int_t | lda, | ||
magmaDoubleComplex * | B, | ||
magma_int_t | ldb, | ||
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 ) |
ZHEGVD computes all the eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x.
Here A and B are assumed to be Hermitian and B is also positive definite. 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.
[in] | ngpu | INTEGER Number of GPUs to use. ngpu > 0. |
[in] | itype | INTEGER Specifies the problem type to be solved:
|
[in] | jobz | magma_vec_t
|
[in] | range | magma_range_t
|
[in] | uplo | magma_uplo_t
|
[in] | n | INTEGER The order of the matrices A and B. N >= 0. |
[in,out] | A | COMPLEX_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 matrix Z of eigenvectors. The eigenvectors are normalized as follows: if ITYPE = 1 or 2, Z**H*B*Z = I; if ITYPE = 3, Z**H*inv(B)*Z = I. If JOBZ = MagmaNoVec, then on exit the upper triangle (if UPLO=MagmaUpper) or the lower triangle (if UPLO=MagmaLower) of A, including the diagonal, is destroyed. |
[in] | lda | INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | COMPLEX_16 array, dimension (LDB, N) On entry, the Hermitian matrix B. If UPLO = MagmaUpper, the leading N-by-N upper triangular part of B contains the upper triangular part of the matrix B. If UPLO = MagmaLower, the leading N-by-N lower triangular part of B contains the lower triangular part of the matrix B. On exit, if INFO <= N, the part of B containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = U**H*U or B = L*L**H. |
[in] | ldb | INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in] | vl | DOUBLE PRECISION |
[in] | vu | DOUBLE 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] | il | INTEGER |
[in] | iu | INTEGER 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] | m | INTEGER The total number of eigenvalues found. 0 <= M <= N. If RANGE = MagmaRangeAll, M = N, and if RANGE = MagmaRangeI, M = IU-IL+1. |
[out] | w | DOUBLE 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] | lwork | INTEGER The length of the array WORK.
|
[out] | rwork | (workspace) DOUBLE PRECISION array, dimension (MAX(1,LRWORK)) On exit, if INFO = 0, RWORK[0] returns the optimal LRWORK. |
[in] | lrwork | INTEGER The dimension of the array RWORK.
|
[out] | iwork | (workspace) INTEGER array, dimension (MAX(1,LIWORK)) On exit, if INFO = 0, IWORK[0] returns the optimal LIWORK. |
[in] | liwork | INTEGER The dimension of the array IWORK.
|
[out] | info | INTEGER
|
Based on contributions by Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
Modified so that no backsubstitution is performed if ZHEEVD fails to converge (NEIG in old code could be greater than N causing out of bounds reference to A - reported by Ralf Meyer). Also corrected the description of INFO and the test on ITYPE. Sven, 16 Feb 05.