Actual source code: dense.c
1: /*
2: This file contains routines for handling small-size dense problems.
3: All routines are simply wrappers to LAPACK routines. Matrices passed in
4: as arguments are assumed to be square matrices stored in column-major
5: format with a leading dimension equal to the number of rows.
6: */
7: #include slepceps.h
8: #include slepcblaslapack.h
12: /*@
13: EPSDenseNHEP - Solves a dense standard non-Hermitian Eigenvalue Problem.
15: Not Collective
17: Input Parameters:
18: + n - dimension of the eigenproblem
19: - A - pointer to the array containing the matrix values
21: Output Parameters:
22: + w - pointer to the array to store the computed eigenvalues
23: . wi - imaginary part of the eigenvalues (only when using real numbers)
24: . V - pointer to the array to store right eigenvectors
25: - W - pointer to the array to store left eigenvectors
27: Notes:
28: If either V or W are PETSC_NULL then the corresponding eigenvectors are
29: not computed.
31: Matrix A is overwritten.
32:
33: This routine uses LAPACK routines xGEEVX.
35: Level: developer
37: .seealso: EPSDenseGNHEP(), EPSDenseHEP(), EPSDenseGHEP()
38: @*/
39: PetscErrorCode EPSDenseNHEP(int n,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
40: {
41: #if defined(SLEPC_MISSING_LAPACK_GEEVX)
43: SETERRQ(PETSC_ERR_SUP,"GEEVX - Lapack routine is unavailable.");
44: #else
46: PetscReal abnrm,*scale,dummy;
47: PetscScalar *work;
48: int ilo,ihi,lwork = 4*n,info;
49: const char *jobvr,*jobvl;
50: #if defined(PETSC_USE_COMPLEX)
51: PetscReal *rwork;
52: #else
53: int idummy;
54: #endif
57: if (V) jobvr = "V";
58: else jobvr = "N";
59: if (W) jobvl = "V";
60: else jobvl = "N";
61: PetscMalloc(lwork*sizeof(PetscScalar),&work);
62: PetscMalloc(n*sizeof(PetscReal),&scale);
63: #if defined(PETSC_USE_COMPLEX)
64: PetscMalloc(2*n*sizeof(PetscReal),&rwork);
65: LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,rwork,&info,1,1,1,1);
66: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZGEEVX %d",info);
67: PetscFree(rwork);
68: #else
69: LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,wi,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,&idummy,&info,1,1,1,1);
70: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DGEEVX %d",info);
71: #endif
72: PetscFree(work);
73: PetscFree(scale);
74: return(0);
75: #endif
76: }
80: /*@
81: EPSDenseGNHEP - Solves a dense Generalized non-Hermitian Eigenvalue Problem.
83: Not Collective
85: Input Parameters:
86: + n - dimension of the eigenproblem
87: . A - pointer to the array containing the matrix values for A
88: - B - pointer to the array containing the matrix values for B
90: Output Parameters:
91: + w - pointer to the array to store the computed eigenvalues
92: . wi - imaginary part of the eigenvalues (only when using real numbers)
93: . V - pointer to the array to store right eigenvectors
94: - W - pointer to the array to store left eigenvectors
96: Notes:
97: If either V or W are PETSC_NULL then the corresponding eigenvectors are
98: not computed.
100: Matrices A and B are overwritten.
101:
102: This routine uses LAPACK routines xGGEVX.
104: Level: developer
106: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGHEP()
107: @*/
108: PetscErrorCode EPSDenseGNHEP(int n,PetscScalar *A,PetscScalar *B,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
109: {
110: #if defined(SLEPC_MISSING_LAPACK_GGEVX)
112: SETERRQ(PETSC_ERR_SUP,"GGEVX - Lapack routine is unavailable.");
113: #else
115: PetscReal *rscale,*lscale,abnrm,bbnrm,dummy;
116: PetscScalar *alpha,*beta,*work;
117: int i,ilo,ihi,idummy,info;
118: const char *jobvr,*jobvl;
119: #if defined(PETSC_USE_COMPLEX)
120: PetscReal *rwork;
121: int lwork = 2*n;
122: #else
123: PetscReal *alphai;
124: int lwork = 6*n;
125: #endif
128: if (V) jobvr = "V";
129: else jobvr = "N";
130: if (W) jobvl = "V";
131: else jobvl = "N";
132: PetscMalloc(n*sizeof(PetscScalar),&alpha);
133: PetscMalloc(n*sizeof(PetscScalar),&beta);
134: PetscMalloc(n*sizeof(PetscReal),&rscale);
135: PetscMalloc(n*sizeof(PetscReal),&lscale);
136: PetscMalloc(lwork*sizeof(PetscScalar),&work);
137: #if defined(PETSC_USE_COMPLEX)
138: PetscMalloc(6*n*sizeof(PetscReal),&rwork);
139: LAPACKggevx_("B",jobvl,jobvr,"N",&n,A,&n,B,&n,alpha,beta,W,&n,V,&n,&ilo,&ihi,lscale,rscale,&abnrm,&bbnrm,&dummy,&dummy,work,&lwork,rwork,&idummy,&idummy,&info,1,1,1,1);
140: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZGGEVX %d",info);
141: for (i=0;i<n;i++) {
142: w[i] = alpha[i]/beta[i];
143: }
144: PetscFree(rwork);
145: #else
146: PetscMalloc(n*sizeof(PetscReal),&alphai);
147: LAPACKggevx_("B",jobvl,jobvr,"N",&n,A,&n,B,&n,alpha,alphai,beta,W,&n,V,&n,&ilo,&ihi,lscale,rscale,&abnrm,&bbnrm,&dummy,&dummy,work,&lwork,&idummy,&idummy,&info,1,1,1,1);
148: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DGGEVX %d",info);
149: for (i=0;i<n;i++) {
150: w[i] = alpha[i]/beta[i];
151: wi[i] = alphai[i]/beta[i];
152: }
153: PetscFree(alphai);
154: #endif
155: PetscFree(alpha);
156: PetscFree(beta);
157: PetscFree(rscale);
158: PetscFree(lscale);
159: PetscFree(work);
160: return(0);
161: #endif
162: }
166: /*@
167: EPSDenseHEP - Solves a dense standard Hermitian Eigenvalue Problem.
169: Not Collective
171: Input Parameters:
172: + n - dimension of the eigenproblem
173: . A - pointer to the array containing the matrix values
174: - lda - leading dimension of A
176: Output Parameters:
177: + w - pointer to the array to store the computed eigenvalues
178: - V - pointer to the array to store the eigenvectors
180: Notes:
181: If V is PETSC_NULL then the eigenvectors are not computed.
183: Matrix A is overwritten.
184:
185: This routine uses LAPACK routines DSYEVR or ZHEEVR.
187: Level: developer
189: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
190: @*/
191: PetscErrorCode EPSDenseHEP(int n,PetscScalar *A,int lda,PetscReal *w,PetscScalar *V)
192: {
193: #if defined(SLEPC_MISSING_LAPACK_SYEVR) || defined(SLEPC_MISSING_LAPACK_HEEVR)
195: SETERRQ(PETSC_ERR_SUP,"DSYEVR/ZHEEVR - Lapack routine is unavailable.");
196: #else
198: PetscReal abstol = 0.0,vl,vu;
199: PetscScalar *work;
200: int il,iu,m,*isuppz,*iwork,liwork = 10*n,info;
201: const char *jobz;
202: #if defined(PETSC_USE_COMPLEX)
203: PetscReal *rwork;
204: int lwork = 18*n,lrwork = 24*n;
205: #else
206: int lwork = 26*n;
207: #endif
210: if (V) jobz = "V";
211: else jobz = "N";
212: PetscMalloc(2*n*sizeof(int),&isuppz);
213: PetscMalloc(lwork*sizeof(PetscScalar),&work);
214: PetscMalloc(liwork*sizeof(int),&iwork);
215: #if defined(PETSC_USE_COMPLEX)
216: PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
217: LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,rwork,&lrwork,iwork,&liwork,&info,1,1,1);
218: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEEVR %d",info);
219: PetscFree(rwork);
220: #else
221: LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1,1);
222: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYEVR %d",info);
223: #endif
224: PetscFree(isuppz);
225: PetscFree(work);
226: PetscFree(iwork);
227: return(0);
228: #endif
229: }
233: /*@
234: EPSDenseGHEP - Solves a dense Generalized Hermitian Eigenvalue Problem.
236: Not Collective
238: Input Parameters:
239: + n - dimension of the eigenproblem
240: . A - pointer to the array containing the matrix values for A
241: - B - pointer to the array containing the matrix values for B
243: Output Parameters:
244: + w - pointer to the array to store the computed eigenvalues
245: - V - pointer to the array to store the eigenvectors
247: Notes:
248: If V is PETSC_NULL then the eigenvectors are not computed.
250: Matrices A and B are overwritten.
251:
252: This routine uses LAPACK routines DSYGVD or ZHEGVD.
254: Level: developer
256: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseHEP()
257: @*/
258: PetscErrorCode EPSDenseGHEP(int n,PetscScalar *A,PetscScalar *B,PetscReal *w,PetscScalar *V)
259: {
260: #if defined(SLEPC_MISSING_LAPACK_SYGVD) || defined(SLEPC_MISSING_LAPACK_HEGVD)
262: SETERRQ(PETSC_ERR_SUP,"DSYGVD/ZHEGVD - Lapack routine is unavailable.");
263: #else
265: PetscScalar *work;
266: int itype = 1,*iwork,info,
267: liwork = V ? 5*n+3 : 1;
268: const char *jobz;
269: #if defined(PETSC_USE_COMPLEX)
270: PetscReal *rwork;
271: int lwork = V ? n*n+2*n : n+1,
272: lrwork = V ? 2*n*n+5*n+1 : n;
273: #else
274: int lwork = V ? 2*n*n+6*n+1 : 2*n+1;
275: #endif
278: if (V) jobz = "V";
279: else jobz = "N";
280: PetscMalloc(lwork*sizeof(PetscScalar),&work);
281: PetscMalloc(liwork*sizeof(int),&iwork);
282: #if defined(PETSC_USE_COMPLEX)
283: PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
284: LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,rwork,&lrwork,iwork,&liwork,&info,1,1);
285: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
286: PetscFree(rwork);
287: #else
288: LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,iwork,&liwork,&info,1,1);
289: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
290: #endif
291: if (V) {
292: PetscMemcpy(V,A,n*n*sizeof(PetscScalar));
293: }
294: PetscFree(work);
295: PetscFree(iwork);
296: return(0);
297: #endif
298: }
302: /*@
303: EPSDenseHessenberg - Computes the Hessenberg form of a dense matrix.
305: Not Collective
307: Input Parameters:
308: + n - dimension of the matrix
309: . k - first active column
310: - lda - leading dimension of A
312: Input/Output Parameters:
313: + A - on entry, the full matrix; on exit, the upper Hessenberg matrix (H)
314: - Q - on exit, orthogonal matrix of vectors A = Q*H*Q'
316: Notes:
317: Only active columns (from k to n) are computed.
319: Both A and Q are overwritten.
320:
321: This routine uses LAPACK routines xGEHRD and xORGHR/xUNGHR.
323: Level: developer
325: .seealso: EPSDenseSchur(), EPSSortDenseSchur(), EPSDenseTridiagonal()
326: @*/
327: PetscErrorCode EPSDenseHessenberg(int n,int k,PetscScalar *A,int lda,PetscScalar *Q)
328: {
329: #if defined(SLEPC_MISSING_LAPACK_GEHRD) || defined(SLEPC_MISSING_LAPACK_ORGHR) || defined(SLEPC_MISSING_LAPACK_UNGHR)
331: SETERRQ(PETSC_ERR_SUP,"GEHRD,ORGHR/UNGHR - Lapack routines are unavailable.");
332: #else
333: PetscScalar *tau,*work;
335: int i,j,ilo,lwork,info;
338: PetscMalloc(n*sizeof(PetscScalar),&tau);
339: lwork = n;
340: PetscMalloc(lwork*sizeof(PetscScalar),&work);
341: ilo = k+1;
342: LAPACKgehrd_(&n,&ilo,&n,A,&lda,tau,work,&lwork,&info);
343: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGEHRD %d",info);
344: for (j=0;j<n-1;j++) {
345: for (i=j+2;i<n;i++) {
346: Q[i+j*n] = A[i+j*lda];
347: A[i+j*lda] = 0.0;
348: }
349: }
350: LAPACKorghr_(&n,&ilo,&n,Q,&n,tau,work,&lwork,&info);
351: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORGHR %d",info);
352: PetscFree(tau);
353: PetscFree(work);
354: return(0);
355: #endif
356: }
360: /*@
361: EPSDenseSchur - Computes the upper (quasi-)triangular form of a dense
362: upper Hessenberg matrix.
364: Not Collective
366: Input Parameters:
367: + n - dimension of the matrix
368: . k - first active column
369: - ldh - leading dimension of H
371: Input/Output Parameters:
372: + H - on entry, the upper Hessenber matrix; on exit, the upper
373: (quasi-)triangular matrix (T)
374: - Z - on entry, initial transformation matrix; on exit, orthogonal
375: matrix of Schur vectors
377: Output Parameters:
378: + wr - pointer to the array to store the computed eigenvalues
379: - wi - imaginary part of the eigenvalues (only when using real numbers)
381: Notes:
382: This function computes the (real) Schur decomposition of an upper
383: Hessenberg matrix H: H*Z = Z*T, where T is an upper (quasi-)triangular
384: matrix (returned in H), and Z is the orthogonal matrix of Schur vectors.
385: Eigenvalues are extracted from the diagonal blocks of T and returned in
386: wr,wi. Transformations are accumulated in Z so that on entry it can
387: contain the transformation matrix associated to the Hessenberg reduction.
389: Only active columns (from k to n) are computed.
391: Both H and Z are overwritten.
392:
393: This routine uses LAPACK routines xHSEQR.
395: Level: developer
397: .seealso: EPSDenseHessenberg(), EPSSortDenseSchur(), EPSDenseTridiagonal()
398: @*/
399: PetscErrorCode EPSDenseSchur(int n,int k,PetscScalar *H,int ldh,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
400: {
401: #if defined(SLEPC_MISSING_LAPACK_HSEQR)
403: SETERRQ(PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable.");
404: #else
406: int ilo,lwork,info;
407: PetscScalar *work;
408: #if !defined(PETSC_USE_COMPLEX)
409: int j;
410: #endif
411:
413: lwork = n;
414: PetscMalloc(lwork*sizeof(PetscScalar),&work);
415: ilo = k+1;
416: #if !defined(PETSC_USE_COMPLEX)
417: LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,wi,Z,&n,work,&lwork,&info,1,1);
418: for (j=0;j<k;j++) {
419: if (j==n-1 || H[j*ldh+j+1] == 0.0) {
420: /* real eigenvalue */
421: wr[j] = H[j*ldh+j];
422: wi[j] = 0.0;
423: } else {
424: /* complex eigenvalue */
425: wr[j] = H[j*ldh+j];
426: wr[j+1] = H[j*ldh+j];
427: wi[j] = sqrt(PetscAbsReal(H[j*ldh+j+1])) *
428: sqrt(PetscAbsReal(H[(j+1)*ldh+j]));
429: wi[j+1] = -wi[j];
430: j++;
431: }
432: }
433: #else
434: LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,Z,&n,work,&lwork,&info,1,1);
435: #endif
436: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);
438: PetscFree(work);
439: return(0);
440: #endif
441: }
445: /*@
446: EPSSortDenseSchur - Reorders the Schur decomposition computed by
447: EPSDenseSchur().
449: Not Collective
451: Input Parameters:
452: + n - dimension of the matrix
453: . k - first active column
454: . ldt - leading dimension of T
455: - which - eigenvalue sort order
457: Input/Output Parameters:
458: + T - the upper (quasi-)triangular matrix
459: - Z - the orthogonal matrix of Schur vectors
461: Output Parameters:
462: + wr - pointer to the array to store the computed eigenvalues
463: - wi - imaginary part of the eigenvalues (only when using real numbers)
465: Notes:
466: This function reorders the eigenvalues in wr,wi located in positions k
467: to n according to the sort order specified in which. The Schur
468: decomposition Z*T*Z^T, is also reordered by means of rotations so that
469: eigenvalues in the diagonal blocks of T follow the same order.
471: Both T and Z are overwritten.
472:
473: This routine uses LAPACK routines xTREXC.
475: Level: developer
477: .seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
478: @*/
479: PetscErrorCode EPSSortDenseSchur(int n,int k,PetscScalar *T,int ldt,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi,EPSWhich which)
480: {
481: #if defined(SLEPC_MISSING_LAPACK_TREXC)
483: SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
484: #else
485: int i,j,ifst,ilst,info,pos;
486: #if !defined(PETSC_USE_COMPLEX)
487: PetscScalar *work;
489: #endif
490: PetscReal value,v;
491:
493: #if !defined(PETSC_USE_COMPLEX)
494: PetscMalloc(n*sizeof(PetscScalar),&work);
495: #endif
496:
497: for (i=k;i<n-1;i++) {
498: switch(which) {
499: case EPS_LARGEST_MAGNITUDE:
500: case EPS_SMALLEST_MAGNITUDE:
501: value = SlepcAbsEigenvalue(wr[i],wi[i]);
502: break;
503: case EPS_LARGEST_REAL:
504: case EPS_SMALLEST_REAL:
505: value = PetscRealPart(wr[i]);
506: break;
507: case EPS_LARGEST_IMAGINARY:
508: case EPS_SMALLEST_IMAGINARY:
509: #if !defined(PETSC_USE_COMPLEX)
510: value = PetscAbsReal(wi[i]);
511: #else
512: value = PetscImaginaryPart(wr[i]);
513: #endif
514: break;
515: default: SETERRQ(1,"Wrong value of which");
516: }
517: pos = 0;
518: for (j=i+1;j<n;j++) {
519: switch(which) {
520: case EPS_LARGEST_MAGNITUDE:
521: case EPS_SMALLEST_MAGNITUDE:
522: v = SlepcAbsEigenvalue(wr[j],wi[j]);
523: break;
524: case EPS_LARGEST_REAL:
525: case EPS_SMALLEST_REAL:
526: v = PetscRealPart(wr[j]);
527: break;
528: case EPS_LARGEST_IMAGINARY:
529: case EPS_SMALLEST_IMAGINARY:
530: #if !defined(PETSC_USE_COMPLEX)
531: v = PetscAbsReal(wi[j]);
532: #else
533: v = PetscImaginaryPart(wr[j]);
534: #endif
535: break;
536: default: SETERRQ(1,"Wrong value of which");
537: }
538: switch(which) {
539: case EPS_LARGEST_MAGNITUDE:
540: case EPS_LARGEST_REAL:
541: case EPS_LARGEST_IMAGINARY:
542: if (v > value) {
543: value = v;
544: pos = j;
545: }
546: break;
547: case EPS_SMALLEST_MAGNITUDE:
548: case EPS_SMALLEST_REAL:
549: case EPS_SMALLEST_IMAGINARY:
550: if (v < value) {
551: value = v;
552: pos = j;
553: }
554: break;
555: default: SETERRQ(1,"Wrong value of which");
556: }
557: #if !defined(PETSC_USE_COMPLEX)
558: if (wi[j] != 0) j++;
559: #endif
560: }
561: if (pos) {
562: ifst = pos + 1;
563: ilst = i + 1;
564: #if !defined(PETSC_USE_COMPLEX)
565: LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,work,&info,1);
566: #else
567: LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,&info,1);
568: #endif
569: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
570:
571: for (j=k;j<n;j++) {
572: #if !defined(PETSC_USE_COMPLEX)
573: if (j==n-1 || T[j*ldt+j+1] == 0.0) {
574: /* real eigenvalue */
575: wr[j] = T[j*ldt+j];
576: wi[j] = 0.0;
577: } else {
578: /* complex eigenvalue */
579: wr[j] = T[j*ldt+j];
580: wr[j+1] = T[j*ldt+j];
581: wi[j] = sqrt(PetscAbsReal(T[j*ldt+j+1])) *
582: sqrt(PetscAbsReal(T[(j+1)*ldt+j]));
583: wi[j+1] = -wi[j];
584: j++;
585: }
586: #else
587: wr[j] = T[j*(ldt+1)];
588: #endif
589: }
590: }
591: #if !defined(PETSC_USE_COMPLEX)
592: if (wi[i] != 0) i++;
593: #endif
594: }
595:
596: #if !defined(PETSC_USE_COMPLEX)
597: PetscFree(work);
598: #endif
599: return(0);
601: #endif
602: }
606: /*@
607: EPSDenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
609: Not Collective
611: Input Parameters:
612: + n - dimension of the eigenproblem
613: . A - pointer to the array containing the matrix values
614: - lda - leading dimension of A
616: Output Parameters:
617: + w - pointer to the array to store the computed eigenvalues
618: - V - pointer to the array to store the eigenvectors
620: Notes:
621: If V is PETSC_NULL then the eigenvectors are not computed.
623: This routine use LAPACK routines DSTEVR.
625: Level: developer
627: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
628: @*/
629: PetscErrorCode EPSDenseTridiagonal(int n,PetscScalar *A,int lda,PetscReal *w,PetscScalar *V)
630: {
631: #if defined(SLEPC_MISSING_LAPACK_STEVR)
633: SETERRQ(PETSC_ERR_SUP,"DSTEVR - Lapack routine is unavailable.");
634: #else
636: PetscReal abstol = 0.0,vl,vu,*D,*E,*work;
637: int i,il,iu,m,*isuppz,lwork = 20*n,*iwork,liwork = 10*n,info;
638: const char *jobz;
639: #if defined(PETSC_USE_COMPLEX)
640: int j;
641: PetscReal *VV;
642: #endif
643:
645: if (V) {
646: jobz = "V";
647: #if defined(PETSC_USE_COMPLEX)
648: PetscMalloc(n*n*sizeof(PetscReal),&VV);
649: #endif
650: } else jobz = "N";
651: PetscMalloc(n*sizeof(PetscReal),&D);
652: PetscMalloc(n*sizeof(PetscReal),&E);
653: PetscMalloc(2*n*sizeof(int),&isuppz);
654: PetscMalloc(lwork*sizeof(PetscReal),&work);
655: PetscMalloc(liwork*sizeof(int),&iwork);
656: for (i=0;i<n-1;i++) {
657: D[i] = PetscRealPart(A[i*(lda+1)]);
658: E[i] = PetscRealPart(A[i*(lda+1)+1]);
659: }
660: D[n-1] = PetscRealPart(A[(n-1)*(lda+1)]);
661: #if defined(PETSC_USE_COMPLEX)
662: LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1);
663: #else
664: LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1);
665: #endif
666: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
667: #if defined(PETSC_USE_COMPLEX)
668: if (V) {
669: for (i=0;i<n;i++)
670: for (j=0;j<n;j++)
671: V[i*n+j] = VV[i*n+j];
672: PetscFree(VV);
673: }
674: #endif
675: PetscFree(D);
676: PetscFree(E);
677: PetscFree(isuppz);
678: PetscFree(work);
679: PetscFree(iwork);
680: return(0);
681: #endif
682: }