00001 /* ssygvd.f -- translated by f2c (version 20061008). 00002 You must link the resulting object file with libf2c: 00003 on Microsoft Windows system, link with libf2c.lib; 00004 on Linux or Unix systems, link with .../path/to/libf2c.a -lm 00005 or, if you install libf2c.a in a standard place, with -lf2c -lm 00006 -- in that order, at the end of the command line, as in 00007 cc *.o -lf2c -lm 00008 Source for libf2c is in /netlib/f2c/libf2c.zip, e.g., 00009 00010 http://www.netlib.org/f2c/libf2c.zip 00011 */ 00012 00013 #include "f2c.h" 00014 #include "blaswrap.h" 00015 00016 /* Table of constant values */ 00017 00018 static real c_b11 = 1.f; 00019 00020 /* Subroutine */ int ssygvd_(integer *itype, char *jobz, char *uplo, integer * 00021 n, real *a, integer *lda, real *b, integer *ldb, real *w, real *work, 00022 integer *lwork, integer *iwork, integer *liwork, integer *info) 00023 { 00024 /* System generated locals */ 00025 integer a_dim1, a_offset, b_dim1, b_offset, i__1; 00026 real r__1, r__2; 00027 00028 /* Local variables */ 00029 integer lopt; 00030 extern logical lsame_(char *, char *); 00031 integer lwmin; 00032 char trans[1]; 00033 integer liopt; 00034 logical upper; 00035 extern /* Subroutine */ int strmm_(char *, char *, char *, char *, 00036 integer *, integer *, real *, real *, integer *, real *, integer * 00037 ); 00038 logical wantz; 00039 extern /* Subroutine */ int strsm_(char *, char *, char *, char *, 00040 integer *, integer *, real *, real *, integer *, real *, integer * 00041 ), xerbla_(char *, integer *); 00042 integer liwmin; 00043 extern /* Subroutine */ int spotrf_(char *, integer *, real *, integer *, 00044 integer *), ssyevd_(char *, char *, integer *, real *, 00045 integer *, real *, real *, integer *, integer *, integer *, 00046 integer *); 00047 logical lquery; 00048 extern /* Subroutine */ int ssygst_(integer *, char *, integer *, real *, 00049 integer *, real *, integer *, integer *); 00050 00051 00052 /* -- LAPACK driver routine (version 3.2) -- */ 00053 /* Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd.. */ 00054 /* November 2006 */ 00055 00056 /* .. Scalar Arguments .. */ 00057 /* .. */ 00058 /* .. Array Arguments .. */ 00059 /* .. */ 00060 00061 /* Purpose */ 00062 /* ======= */ 00063 00064 /* SSYGVD computes all the eigenvalues, and optionally, the eigenvectors */ 00065 /* of a real generalized symmetric-definite eigenproblem, of the form */ 00066 /* A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. Here A and */ 00067 /* B are assumed to be symmetric and B is also positive definite. */ 00068 /* If eigenvectors are desired, it uses a divide and conquer algorithm. */ 00069 00070 /* The divide and conquer algorithm makes very mild assumptions about */ 00071 /* floating point arithmetic. It will work on machines with a guard */ 00072 /* digit in add/subtract, or on those binary machines without guard */ 00073 /* digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or */ 00074 /* Cray-2. It could conceivably fail on hexadecimal or decimal machines */ 00075 /* without guard digits, but we know of none. */ 00076 00077 /* Arguments */ 00078 /* ========= */ 00079 00080 /* ITYPE (input) INTEGER */ 00081 /* Specifies the problem type to be solved: */ 00082 /* = 1: A*x = (lambda)*B*x */ 00083 /* = 2: A*B*x = (lambda)*x */ 00084 /* = 3: B*A*x = (lambda)*x */ 00085 00086 /* JOBZ (input) CHARACTER*1 */ 00087 /* = 'N': Compute eigenvalues only; */ 00088 /* = 'V': Compute eigenvalues and eigenvectors. */ 00089 00090 /* UPLO (input) CHARACTER*1 */ 00091 /* = 'U': Upper triangles of A and B are stored; */ 00092 /* = 'L': Lower triangles of A and B are stored. */ 00093 00094 /* N (input) INTEGER */ 00095 /* The order of the matrices A and B. N >= 0. */ 00096 00097 /* A (input/output) REAL array, dimension (LDA, N) */ 00098 /* On entry, the symmetric matrix A. If UPLO = 'U', the */ 00099 /* leading N-by-N upper triangular part of A contains the */ 00100 /* upper triangular part of the matrix A. If UPLO = 'L', */ 00101 /* the leading N-by-N lower triangular part of A contains */ 00102 /* the lower triangular part of the matrix A. */ 00103 00104 /* On exit, if JOBZ = 'V', then if INFO = 0, A contains the */ 00105 /* matrix Z of eigenvectors. The eigenvectors are normalized */ 00106 /* as follows: */ 00107 /* if ITYPE = 1 or 2, Z**T*B*Z = I; */ 00108 /* if ITYPE = 3, Z**T*inv(B)*Z = I. */ 00109 /* If JOBZ = 'N', then on exit the upper triangle (if UPLO='U') */ 00110 /* or the lower triangle (if UPLO='L') of A, including the */ 00111 /* diagonal, is destroyed. */ 00112 00113 /* LDA (input) INTEGER */ 00114 /* The leading dimension of the array A. LDA >= max(1,N). */ 00115 00116 /* B (input/output) REAL array, dimension (LDB, N) */ 00117 /* On entry, the symmetric matrix B. If UPLO = 'U', the */ 00118 /* leading N-by-N upper triangular part of B contains the */ 00119 /* upper triangular part of the matrix B. If UPLO = 'L', */ 00120 /* the leading N-by-N lower triangular part of B contains */ 00121 /* the lower triangular part of the matrix B. */ 00122 00123 /* On exit, if INFO <= N, the part of B containing the matrix is */ 00124 /* overwritten by the triangular factor U or L from the Cholesky */ 00125 /* factorization B = U**T*U or B = L*L**T. */ 00126 00127 /* LDB (input) INTEGER */ 00128 /* The leading dimension of the array B. LDB >= max(1,N). */ 00129 00130 /* W (output) REAL array, dimension (N) */ 00131 /* If INFO = 0, the eigenvalues in ascending order. */ 00132 00133 /* WORK (workspace/output) REAL array, dimension (MAX(1,LWORK)) */ 00134 /* On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */ 00135 00136 /* LWORK (input) INTEGER */ 00137 /* The dimension of the array WORK. */ 00138 /* If N <= 1, LWORK >= 1. */ 00139 /* If JOBZ = 'N' and N > 1, LWORK >= 2*N+1. */ 00140 /* If JOBZ = 'V' and N > 1, LWORK >= 1 + 6*N + 2*N**2. */ 00141 00142 /* If LWORK = -1, then a workspace query is assumed; the routine */ 00143 /* only calculates the optimal sizes of the WORK and IWORK */ 00144 /* arrays, returns these values as the first entries of the WORK */ 00145 /* and IWORK arrays, and no error message related to LWORK or */ 00146 /* LIWORK is issued by XERBLA. */ 00147 00148 /* IWORK (workspace/output) INTEGER array, dimension (MAX(1,LIWORK)) */ 00149 /* On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK. */ 00150 00151 /* LIWORK (input) INTEGER */ 00152 /* The dimension of the array IWORK. */ 00153 /* If N <= 1, LIWORK >= 1. */ 00154 /* If JOBZ = 'N' and N > 1, LIWORK >= 1. */ 00155 /* If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N. */ 00156 00157 /* If LIWORK = -1, then a workspace query is assumed; the */ 00158 /* routine only calculates the optimal sizes of the WORK and */ 00159 /* IWORK arrays, returns these values as the first entries of */ 00160 /* the WORK and IWORK arrays, and no error message related to */ 00161 /* LWORK or LIWORK is issued by XERBLA. */ 00162 00163 /* INFO (output) INTEGER */ 00164 /* = 0: successful exit */ 00165 /* < 0: if INFO = -i, the i-th argument had an illegal value */ 00166 /* > 0: SPOTRF or SSYEVD returned an error code: */ 00167 /* <= N: if INFO = i and JOBZ = 'N', then the algorithm */ 00168 /* failed to converge; i off-diagonal elements of an */ 00169 /* intermediate tridiagonal form did not converge to */ 00170 /* zero; */ 00171 /* if INFO = i and JOBZ = 'V', then the algorithm */ 00172 /* failed to compute an eigenvalue while working on */ 00173 /* the submatrix lying in rows and columns INFO/(N+1) */ 00174 /* through mod(INFO,N+1); */ 00175 /* > N: if INFO = N + i, for 1 <= i <= N, then the leading */ 00176 /* minor of order i of B is not positive definite. */ 00177 /* The factorization of B could not be completed and */ 00178 /* no eigenvalues or eigenvectors were computed. */ 00179 00180 /* Further Details */ 00181 /* =============== */ 00182 00183 /* Based on contributions by */ 00184 /* Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA */ 00185 00186 /* Modified so that no backsubstitution is performed if SSYEVD fails to */ 00187 /* converge (NEIG in old code could be greater than N causing out of */ 00188 /* bounds reference to A - reported by Ralf Meyer). Also corrected the */ 00189 /* description of INFO and the test on ITYPE. Sven, 16 Feb 05. */ 00190 /* ===================================================================== */ 00191 00192 /* .. Parameters .. */ 00193 /* .. */ 00194 /* .. Local Scalars .. */ 00195 /* .. */ 00196 /* .. External Functions .. */ 00197 /* .. */ 00198 /* .. External Subroutines .. */ 00199 /* .. */ 00200 /* .. Intrinsic Functions .. */ 00201 /* .. */ 00202 /* .. Executable Statements .. */ 00203 00204 /* Test the input parameters. */ 00205 00206 /* Parameter adjustments */ 00207 a_dim1 = *lda; 00208 a_offset = 1 + a_dim1; 00209 a -= a_offset; 00210 b_dim1 = *ldb; 00211 b_offset = 1 + b_dim1; 00212 b -= b_offset; 00213 --w; 00214 --work; 00215 --iwork; 00216 00217 /* Function Body */ 00218 wantz = lsame_(jobz, "V"); 00219 upper = lsame_(uplo, "U"); 00220 lquery = *lwork == -1 || *liwork == -1; 00221 00222 *info = 0; 00223 if (*n <= 1) { 00224 liwmin = 1; 00225 lwmin = 1; 00226 } else if (wantz) { 00227 liwmin = *n * 5 + 3; 00228 /* Computing 2nd power */ 00229 i__1 = *n; 00230 lwmin = *n * 6 + 1 + (i__1 * i__1 << 1); 00231 } else { 00232 liwmin = 1; 00233 lwmin = (*n << 1) + 1; 00234 } 00235 lopt = lwmin; 00236 liopt = liwmin; 00237 if (*itype < 1 || *itype > 3) { 00238 *info = -1; 00239 } else if (! (wantz || lsame_(jobz, "N"))) { 00240 *info = -2; 00241 } else if (! (upper || lsame_(uplo, "L"))) { 00242 *info = -3; 00243 } else if (*n < 0) { 00244 *info = -4; 00245 } else if (*lda < max(1,*n)) { 00246 *info = -6; 00247 } else if (*ldb < max(1,*n)) { 00248 *info = -8; 00249 } 00250 00251 if (*info == 0) { 00252 work[1] = (real) lopt; 00253 iwork[1] = liopt; 00254 00255 if (*lwork < lwmin && ! lquery) { 00256 *info = -11; 00257 } else if (*liwork < liwmin && ! lquery) { 00258 *info = -13; 00259 } 00260 } 00261 00262 if (*info != 0) { 00263 i__1 = -(*info); 00264 xerbla_("SSYGVD", &i__1); 00265 return 0; 00266 } else if (lquery) { 00267 return 0; 00268 } 00269 00270 /* Quick return if possible */ 00271 00272 if (*n == 0) { 00273 return 0; 00274 } 00275 00276 /* Form a Cholesky factorization of B. */ 00277 00278 spotrf_(uplo, n, &b[b_offset], ldb, info); 00279 if (*info != 0) { 00280 *info = *n + *info; 00281 return 0; 00282 } 00283 00284 /* Transform problem to standard eigenvalue problem and solve. */ 00285 00286 ssygst_(itype, uplo, n, &a[a_offset], lda, &b[b_offset], ldb, info); 00287 ssyevd_(jobz, uplo, n, &a[a_offset], lda, &w[1], &work[1], lwork, &iwork[ 00288 1], liwork, info); 00289 /* Computing MAX */ 00290 r__1 = (real) lopt; 00291 lopt = dmax(r__1,work[1]); 00292 /* Computing MAX */ 00293 r__1 = (real) liopt, r__2 = (real) iwork[1]; 00294 liopt = dmax(r__1,r__2); 00295 00296 if (wantz && *info == 0) { 00297 00298 /* Backtransform eigenvectors to the original problem. */ 00299 00300 if (*itype == 1 || *itype == 2) { 00301 00302 /* For A*x=(lambda)*B*x and A*B*x=(lambda)*x; */ 00303 /* backtransform eigenvectors: x = inv(L)'*y or inv(U)*y */ 00304 00305 if (upper) { 00306 *(unsigned char *)trans = 'N'; 00307 } else { 00308 *(unsigned char *)trans = 'T'; 00309 } 00310 00311 strsm_("Left", uplo, trans, "Non-unit", n, n, &c_b11, &b[b_offset] 00312 , ldb, &a[a_offset], lda); 00313 00314 } else if (*itype == 3) { 00315 00316 /* For B*A*x=(lambda)*x; */ 00317 /* backtransform eigenvectors: x = L*y or U'*y */ 00318 00319 if (upper) { 00320 *(unsigned char *)trans = 'T'; 00321 } else { 00322 *(unsigned char *)trans = 'N'; 00323 } 00324 00325 strmm_("Left", uplo, trans, "Non-unit", n, n, &c_b11, &b[b_offset] 00326 , ldb, &a[a_offset], lda); 00327 } 00328 } 00329 00330 work[1] = (real) lopt; 00331 iwork[1] = liopt; 00332 00333 return 0; 00334 00335 /* End of SSYGVD */ 00336 00337 } /* ssygvd_ */