table of contents
dsyevd_2stage.f(3) | LAPACK | dsyevd_2stage.f(3) |
NAME¶
dsyevd_2stage.f
SYNOPSIS¶
Functions/Subroutines¶
subroutine dsyevd_2stage (JOBZ, UPLO, N, A,
LDA, W, WORK, LWORK, IWORK, LIWORK, INFO)
DSYEVD_2STAGE computes the eigenvalues and, optionally, the left and/or
right eigenvectors for SY matrices
Function/Subroutine Documentation¶
subroutine dsyevd_2stage (character JOBZ, character UPLO, integer N, double precision, dimension( lda, * ) A, integer LDA, double precision, dimension( * ) W, double precision, dimension( * ) WORK, integer LWORK, integer, dimension( * ) IWORK, integer LIWORK, integer INFO)¶
DSYEVD_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices
Purpose:
DSYEVD_2STAGE computes all eigenvalues and, optionally, eigenvectors of a
real symmetric matrix A using the 2stage technique for
the reduction to tridiagonal. 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:
JOBZ
JOBZ is CHARACTER*1
= 'N': Compute eigenvalues only;
= 'V': Compute eigenvalues and eigenvectors.
Not available in this release.
UPLO
UPLO is CHARACTER*1
= 'U': Upper triangle of A is stored;
= 'L': Lower triangle of A is stored.
N
N is INTEGER
The order of the matrix A. N >= 0.
A
A is DOUBLE PRECISION array, dimension (LDA, N)
On entry, the symmetric matrix A. If UPLO = 'U', the
leading N-by-N upper triangular part of A contains the
upper triangular part of the matrix A. If UPLO = 'L',
the leading N-by-N lower triangular part of A contains
the lower triangular part of the matrix A.
On exit, if JOBZ = 'V', then if INFO = 0, A contains the
orthonormal eigenvectors of the matrix A.
If JOBZ = 'N', then on exit the lower triangle (if UPLO='L')
or the upper triangle (if UPLO='U') of A, including the
diagonal, is destroyed.
LDA
LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).
W
W is DOUBLE PRECISION array, dimension (N)
If INFO = 0, the eigenvalues in ascending order.
WORK
WORK is DOUBLE PRECISION array,
dimension (LWORK)
On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
LWORK
LWORK is INTEGER
The dimension of the array WORK.
If N <= 1, LWORK must be at least 1.
If JOBZ = 'N' and N > 1, LWORK must be queried.
LWORK = MAX(1, dimension) where
dimension = max(stage1,stage2) + (KD+1)*N + 2*N+1
= N*KD + N*max(KD+1,FACTOPTNB)
+ max(2*KD*KD, KD*NTHREADS)
+ (KD+1)*N + 2*N+1
where KD is the blocking size of the reduction,
FACTOPTNB is the blocking used by the QR or LQ
algorithm, usually FACTOPTNB=128 is a good choice
NTHREADS is the number of threads used when
openMP compilation is enabled, otherwise =1.
If JOBZ = 'V' and N > 1, LWORK must be at least
1 + 6*N + 2*N**2.
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.
IWORK
IWORK is INTEGER array, dimension (MAX(1,LIWORK))
On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
LIWORK
LIWORK is INTEGER
The dimension of the array IWORK.
If N <= 1, LIWORK must be at least 1.
If JOBZ = 'N' and N > 1, LIWORK must be at least 1.
If JOBZ = 'V' and N > 1, LIWORK must be at least 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.
INFO
INFO is INTEGER
= 0: successful exit
< 0: if INFO = -i, the i-th argument had an illegal value
> 0: if INFO = i and JOBZ = 'N', 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 = 'V', 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).
Author:
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date:
November 2017
Contributors:
Jeff Rutter, Computer Science Division, University of
California at Berkeley, USA
Modified by Francoise Tisseur, University of Tennessee
Modified description of INFO. Sven, 16 Feb 05.
Modified by Francoise Tisseur, University of Tennessee
Modified description of INFO. Sven, 16 Feb 05.
Further Details:
All details about the 2stage techniques are available in:
Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
Parallel reduction to condensed forms for symmetric eigenvalue problems
using aggregated fine-grained and memory-aware kernels. In Proceedings
of 2011 International Conference for High Performance Computing,
Networking, Storage and Analysis (SC '11), New York, NY, USA,
Article 8 , 11 pages.
http://doi.acm.org/10.1145/2063384.2063394
A. Haidar, J. Kurzak, P. Luszczek, 2013.
An improved parallel singular value algorithm and its implementation
for multicore hardware, In Proceedings of 2013 International Conference
for High Performance Computing, Networking, Storage and Analysis (SC '13).
Denver, Colorado, USA, 2013.
Article 90, 12 pages.
http://doi.acm.org/10.1145/2503210.2503292
A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
A novel hybrid CPU-GPU generalized eigensolver for electronic structure
calculations based on fine-grained memory aware tasks.
International Journal of High Performance Computing Applications.
Volume 28 Issue 2, Pages 196-209, May 2014.
http://hpc.sagepub.com/content/28/2/196
Definition at line 229 of file dsyevd_2stage.f.
Author¶
Generated automatically by Doxygen for LAPACK from the source code.
Tue Nov 14 2017 | Version 3.8.0 |