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.
7: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
8: SLEPc - Scalable Library for Eigenvalue Problem Computations
9: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
11: This file is part of SLEPc.
12:
13: SLEPc is free software: you can redistribute it and/or modify it under the
14: terms of version 3 of the GNU Lesser General Public License as published by
15: the Free Software Foundation.
17: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
18: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
20: more details.
22: You should have received a copy of the GNU Lesser General Public License
23: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
24: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25: */
27: #include <private/epsimpl.h> /*I "slepceps.h" I*/
28: #include <slepcblaslapack.h>
32: /*@
33: EPSDenseNHEP - Solves a dense standard non-Hermitian Eigenvalue Problem.
35: Not Collective
37: Input Parameters:
38: + n - dimension of the eigenproblem
39: - A - pointer to the array containing the matrix values
41: Output Parameters:
42: + w - pointer to the array to store the computed eigenvalues
43: . wi - imaginary part of the eigenvalues (only when using real numbers)
44: . V - pointer to the array to store right eigenvectors
45: - W - pointer to the array to store left eigenvectors
47: Notes:
48: If either V or W are PETSC_NULL then the corresponding eigenvectors are
49: not computed.
51: Matrix A is overwritten.
52:
53: This routine uses LAPACK routines xGEEVX.
55: Level: developer
57: .seealso: EPSDenseGNHEP(), EPSDenseHEP(), EPSDenseGHEP()
58: @*/
59: PetscErrorCode EPSDenseNHEP(PetscInt n_,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
60: {
61: #if defined(SLEPC_MISSING_LAPACK_GEEVX)
63: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEEVX - Lapack routine is unavailable.");
64: #else
66: PetscReal abnrm,*scale,dummy;
67: PetscScalar *work;
68: PetscBLASInt ilo,ihi,n,lwork,info;
69: const char *jobvr,*jobvl;
70: #if defined(PETSC_USE_COMPLEX)
71: PetscReal *rwork;
72: #else
73: PetscBLASInt idummy;
74: #endif
77: PetscLogEventBegin(EPS_Dense,0,0,0,0);
78: n = PetscBLASIntCast(n_);
79: lwork = PetscBLASIntCast(4*n_);
80: if (V) jobvr = "V";
81: else jobvr = "N";
82: if (W) jobvl = "V";
83: else jobvl = "N";
84: PetscMalloc(lwork*sizeof(PetscScalar),&work);
85: PetscMalloc(n*sizeof(PetscReal),&scale);
86: #if defined(PETSC_USE_COMPLEX)
87: PetscMalloc(2*n*sizeof(PetscReal),&rwork);
88: LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,rwork,&info);
89: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZGEEVX %d",info);
90: PetscFree(rwork);
91: #else
92: LAPACKgeevx_("B",jobvl,jobvr,"N",&n,A,&n,w,wi,W,&n,V,&n,&ilo,&ihi,scale,&abnrm,&dummy,&dummy,work,&lwork,&idummy,&info);
93: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DGEEVX %d",info);
94: #endif
95: PetscFree(work);
96: PetscFree(scale);
97: PetscLogEventEnd(EPS_Dense,0,0,0,0);
98: return(0);
99: #endif
100: }
104: /*@
105: EPSDenseGNHEP - Solves a dense Generalized non-Hermitian Eigenvalue Problem.
107: Not Collective
109: Input Parameters:
110: + n - dimension of the eigenproblem
111: . A - pointer to the array containing the matrix values for A
112: - B - pointer to the array containing the matrix values for B
114: Output Parameters:
115: + w - pointer to the array to store the computed eigenvalues
116: . wi - imaginary part of the eigenvalues (only when using real numbers)
117: . V - pointer to the array to store right eigenvectors
118: - W - pointer to the array to store left eigenvectors
120: Notes:
121: If either V or W are PETSC_NULL then the corresponding eigenvectors are
122: not computed.
124: Matrices A and B are overwritten.
125:
126: This routine uses LAPACK routines xGGEVX.
128: Level: developer
130: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGHEP()
131: @*/
132: PetscErrorCode EPSDenseGNHEP(PetscInt n_,PetscScalar *A,PetscScalar *B,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
133: {
134: #if defined(SLEPC_MISSING_LAPACK_GGEVX)
136: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GGEVX - Lapack routine is unavailable.");
137: #else
139: PetscReal *rscale,*lscale,abnrm,bbnrm,dummy;
140: PetscScalar *alpha,*beta,*work;
141: PetscInt i;
142: PetscBLASInt ilo,ihi,idummy,info,n;
143: const char *jobvr,*jobvl;
144: #if defined(PETSC_USE_COMPLEX)
145: PetscReal *rwork;
146: PetscBLASInt lwork;
147: #else
148: PetscReal *alphai;
149: PetscBLASInt lwork;
150: #endif
153: PetscLogEventBegin(EPS_Dense,0,0,0,0);
154: n = PetscBLASIntCast(n_);
155: #if defined(PETSC_USE_COMPLEX)
156: lwork = PetscBLASIntCast(2*n_);
157: #else
158: lwork = PetscBLASIntCast(6*n_);
159: #endif
160: if (V) jobvr = "V";
161: else jobvr = "N";
162: if (W) jobvl = "V";
163: else jobvl = "N";
164: PetscMalloc(n*sizeof(PetscScalar),&alpha);
165: PetscMalloc(n*sizeof(PetscScalar),&beta);
166: PetscMalloc(n*sizeof(PetscReal),&rscale);
167: PetscMalloc(n*sizeof(PetscReal),&lscale);
168: PetscMalloc(lwork*sizeof(PetscScalar),&work);
169: #if defined(PETSC_USE_COMPLEX)
170: PetscMalloc(6*n*sizeof(PetscReal),&rwork);
171: 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);
172: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZGGEVX %d",info);
173: for (i=0;i<n;i++) {
174: w[i] = alpha[i]/beta[i];
175: }
176: PetscFree(rwork);
177: #else
178: PetscMalloc(n*sizeof(PetscReal),&alphai);
179: 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);
180: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DGGEVX %d",info);
181: for (i=0;i<n;i++) {
182: w[i] = alpha[i]/beta[i];
183: wi[i] = alphai[i]/beta[i];
184: }
185: PetscFree(alphai);
186: #endif
187: PetscFree(alpha);
188: PetscFree(beta);
189: PetscFree(rscale);
190: PetscFree(lscale);
191: PetscFree(work);
192: PetscLogEventEnd(EPS_Dense,0,0,0,0);
193: return(0);
194: #endif
195: }
199: /*@
200: EPSDenseHEP - Solves a dense standard Hermitian Eigenvalue Problem.
202: Not Collective
204: Input Parameters:
205: + n - dimension of the eigenproblem
206: . A - pointer to the array containing the matrix values
207: - lda - leading dimension of A
209: Output Parameters:
210: + w - pointer to the array to store the computed eigenvalues
211: - V - pointer to the array to store the eigenvectors
213: Notes:
214: If V is PETSC_NULL then the eigenvectors are not computed.
216: Matrix A is overwritten.
217:
218: This routine uses LAPACK routines DSYEVR or ZHEEVR.
220: Level: developer
222: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
223: @*/
224: PetscErrorCode EPSDenseHEP(PetscInt n_,PetscScalar *A,PetscInt lda_,PetscReal *w,PetscScalar *V)
225: {
226: #if defined(SLEPC_MISSING_LAPACK_SYEVR) || defined(SLEPC_MISSING_LAPACK_HEEVR)
228: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSYEVR/ZHEEVR - Lapack routine is unavailable.");
229: #else
231: PetscReal abstol = 0.0,vl,vu;
232: PetscScalar *work;
233: PetscBLASInt il,iu,m,*isuppz,*iwork,n,lda,liwork,info;
234: const char *jobz;
235: #if defined(PETSC_USE_COMPLEX)
236: PetscReal *rwork;
237: PetscBLASInt lwork,lrwork;
238: #else
239: PetscBLASInt lwork;
240: #endif
243: PetscLogEventBegin(EPS_Dense,0,0,0,0);
244: n = PetscBLASIntCast(n_);
245: lda = PetscBLASIntCast(lda_);
246: liwork = PetscBLASIntCast(10*n_);
247: #if defined(PETSC_USE_COMPLEX)
248: lwork = PetscBLASIntCast(18*n_);
249: lrwork = PetscBLASIntCast(24*n_);
250: #else
251: lwork = PetscBLASIntCast(26*n_);
252: #endif
253: if (V) jobz = "V";
254: else jobz = "N";
255: PetscMalloc(2*n*sizeof(PetscBLASInt),&isuppz);
256: PetscMalloc(lwork*sizeof(PetscScalar),&work);
257: PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
258: #if defined(PETSC_USE_COMPLEX)
259: PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
260: LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,rwork,&lrwork,iwork,&liwork,&info);
261: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZHEEVR %d",info);
262: PetscFree(rwork);
263: #else
264: LAPACKsyevr_(jobz,"A","L",&n,A,&lda,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info);
265: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSYEVR %d",info);
266: #endif
267: PetscFree(isuppz);
268: PetscFree(work);
269: PetscFree(iwork);
270: PetscLogEventEnd(EPS_Dense,0,0,0,0);
271: return(0);
272: #endif
273: }
277: /*@
278: EPSDenseGHEP - Solves a dense Generalized Hermitian Eigenvalue Problem.
280: Not Collective
282: Input Parameters:
283: + n - dimension of the eigenproblem
284: . A - pointer to the array containing the matrix values for A
285: - B - pointer to the array containing the matrix values for B
287: Output Parameters:
288: + w - pointer to the array to store the computed eigenvalues
289: - V - pointer to the array to store the eigenvectors
291: Notes:
292: If V is PETSC_NULL then the eigenvectors are not computed.
294: Matrices A and B are overwritten.
295:
296: This routine uses LAPACK routines DSYGVD or ZHEGVD.
298: Level: developer
300: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseHEP()
301: @*/
302: PetscErrorCode EPSDenseGHEP(PetscInt n_,PetscScalar *A,PetscScalar *B,PetscReal *w,PetscScalar *V)
303: {
304: #if defined(SLEPC_MISSING_LAPACK_SYGVD) || defined(SLEPC_MISSING_LAPACK_HEGVD)
306: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSYGVD/ZHEGVD - Lapack routine is unavailable.");
307: #else
309: PetscScalar *work;
310: PetscBLASInt itype = 1,*iwork,info,n,liwork;
311: const char *jobz;
312: #if defined(PETSC_USE_COMPLEX)
313: PetscReal *rwork;
314: PetscBLASInt lwork,lrwork;
315: #else
316: PetscBLASInt lwork;
317: #endif
320: PetscLogEventBegin(EPS_Dense,0,0,0,0);
321: n = PetscBLASIntCast(n_);
322: if (V) {
323: jobz = "V";
324: liwork = PetscBLASIntCast(5*n_+3);
325: #if defined(PETSC_USE_COMPLEX)
326: lwork = PetscBLASIntCast(n_*n_+2*n_);
327: lrwork = PetscBLASIntCast(2*n_*n_+5*n_+1);
328: #else
329: lwork = PetscBLASIntCast(2*n_*n_+6*n_+1);
330: #endif
331: } else {
332: jobz = "N";
333: liwork = 1;
334: #if defined(PETSC_USE_COMPLEX)
335: lwork = PetscBLASIntCast(n_+1);
336: lrwork = PetscBLASIntCast(n_);
337: #else
338: lwork = PetscBLASIntCast(2*n_+1);
339: #endif
340: }
341: PetscMalloc(lwork*sizeof(PetscScalar),&work);
342: PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
343: #if defined(PETSC_USE_COMPLEX)
344: PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
345: LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,rwork,&lrwork,iwork,&liwork,&info);
346: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
347: PetscFree(rwork);
348: #else
349: LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,iwork,&liwork,&info);
350: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
351: #endif
352: if (V) {
353: PetscMemcpy(V,A,n*n*sizeof(PetscScalar));
354: }
355: PetscFree(work);
356: PetscFree(iwork);
357: PetscLogEventEnd(EPS_Dense,0,0,0,0);
358: return(0);
359: #endif
360: }
364: /*@
365: EPSDenseHessenberg - Computes the Hessenberg form of a dense matrix.
367: Not Collective
369: Input Parameters:
370: + n - dimension of the matrix
371: . k - first active column
372: - lda - leading dimension of A
374: Input/Output Parameters:
375: + A - on entry, the full matrix; on exit, the upper Hessenberg matrix (H)
376: - Q - on exit, orthogonal matrix of vectors A = Q*H*Q'
378: Notes:
379: Only active columns (from k to n) are computed.
381: Both A and Q are overwritten.
382:
383: This routine uses LAPACK routines xGEHRD and xORGHR/xUNGHR.
385: Level: developer
387: .seealso: EPSDenseSchur(), EPSSortDenseSchur(), EPSDenseTridiagonal()
388: @*/
389: PetscErrorCode EPSDenseHessenberg(PetscInt n_,PetscInt k,PetscScalar *A,PetscInt lda_,PetscScalar *Q)
390: {
391: #if defined(SLEPC_MISSING_LAPACK_GEHRD) || defined(SLEPC_MISSING_LAPACK_ORGHR) || defined(SLEPC_MISSING_LAPACK_UNGHR)
393: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEHRD,ORGHR/UNGHR - Lapack routines are unavailable.");
394: #else
395: PetscScalar *tau,*work;
397: PetscInt i,j;
398: PetscBLASInt ilo,lwork,info,n,lda;
401: PetscLogEventBegin(EPS_Dense,0,0,0,0);
402: n = PetscBLASIntCast(n_);
403: lda = PetscBLASIntCast(lda_);
404: PetscMalloc(n*sizeof(PetscScalar),&tau);
405: lwork = n;
406: PetscMalloc(lwork*sizeof(PetscScalar),&work);
407: ilo = PetscBLASIntCast(k+1);
408: LAPACKgehrd_(&n,&ilo,&n,A,&lda,tau,work,&lwork,&info);
409: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEHRD %d",info);
410: for (j=0;j<n-1;j++) {
411: for (i=j+2;i<n;i++) {
412: Q[i+j*n] = A[i+j*lda];
413: A[i+j*lda] = 0.0;
414: }
415: }
416: LAPACKorghr_(&n,&ilo,&n,Q,&n,tau,work,&lwork,&info);
417: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGHR %d",info);
418: PetscFree(tau);
419: PetscFree(work);
420: PetscLogEventEnd(EPS_Dense,0,0,0,0);
421: return(0);
422: #endif
423: }
427: /*@
428: EPSDenseSchur - Computes the upper (quasi-)triangular form of a dense
429: upper Hessenberg matrix.
431: Not Collective
433: Input Parameters:
434: + n - dimension of the matrix
435: . k - first active column
436: - ldh - leading dimension of H
438: Input/Output Parameters:
439: + H - on entry, the upper Hessenber matrix; on exit, the upper
440: (quasi-)triangular matrix (T)
441: - Z - on entry, initial transformation matrix; on exit, orthogonal
442: matrix of Schur vectors
444: Output Parameters:
445: + wr - pointer to the array to store the computed eigenvalues
446: - wi - imaginary part of the eigenvalues (only when using real numbers)
448: Notes:
449: This function computes the (real) Schur decomposition of an upper
450: Hessenberg matrix H: H*Z = Z*T, where T is an upper (quasi-)triangular
451: matrix (returned in H), and Z is the orthogonal matrix of Schur vectors.
452: Eigenvalues are extracted from the diagonal blocks of T and returned in
453: wr,wi. Transformations are accumulated in Z so that on entry it can
454: contain the transformation matrix associated to the Hessenberg reduction.
456: Only active columns (from k to n) are computed.
458: Both H and Z are overwritten.
459:
460: This routine uses LAPACK routines xHSEQR.
462: Level: developer
464: .seealso: EPSDenseHessenberg(), EPSSortDenseSchur(), EPSDenseTridiagonal()
465: @*/
466: PetscErrorCode EPSDenseSchur(PetscInt n_,PetscInt k,PetscScalar *H,PetscInt ldh_,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
467: {
468: #if defined(SLEPC_MISSING_LAPACK_HSEQR)
470: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable.");
471: #else
473: PetscBLASInt ilo,lwork,info,n,ldh;
474: PetscScalar *work;
475: #if !defined(PETSC_USE_COMPLEX)
476: PetscInt j;
477: #endif
478:
480: PetscLogEventBegin(EPS_Dense,0,0,0,0);
481: n = PetscBLASIntCast(n_);
482: ldh = PetscBLASIntCast(ldh_);
483: lwork = n;
484: PetscMalloc(lwork*sizeof(PetscScalar),&work);
485: ilo = PetscBLASIntCast(k+1);
486: #if !defined(PETSC_USE_COMPLEX)
487: LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,wi,Z,&n,work,&lwork,&info);
488: for (j=0;j<k;j++) {
489: if (j==n-1 || H[j*ldh+j+1] == 0.0) {
490: /* real eigenvalue */
491: wr[j] = H[j*ldh+j];
492: wi[j] = 0.0;
493: } else {
494: /* complex eigenvalue */
495: wr[j] = H[j*ldh+j];
496: wr[j+1] = H[j*ldh+j];
497: wi[j] = PetscSqrtReal(PetscAbsReal(H[j*ldh+j+1])) *
498: PetscSqrtReal(PetscAbsReal(H[(j+1)*ldh+j]));
499: wi[j+1] = -wi[j];
500: j++;
501: }
502: }
503: #else
504: LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,Z,&n,work,&lwork,&info);
505: #endif
506: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);
508: PetscFree(work);
509: PetscLogEventEnd(EPS_Dense,0,0,0,0);
510: return(0);
511: #endif
512: }
516: /*@
517: EPSSortDenseSchur - Reorders the Schur decomposition computed by
518: EPSDenseSchur().
520: Not Collective
522: Input Parameters:
523: + eps - the eigensolver context
524: . n - dimension of the matrix
525: . k - first active column
526: - ldt - leading dimension of T
528: Input/Output Parameters:
529: + T - the upper (quasi-)triangular matrix
530: . Q - the orthogonal matrix of Schur vectors
531: . wr - pointer to the array to store the computed eigenvalues
532: - wi - imaginary part of the eigenvalues (only when using real numbers)
534: Notes:
535: This function reorders the eigenvalues in wr,wi located in positions k
536: to n according to the sort order specified in EPSetWhicheigenpairs.
537: The Schur decomposition Z*T*Z^T, is also reordered by means of rotations
538: so that eigenvalues in the diagonal blocks of T follow the same order.
540: Both T and Q are overwritten.
541:
542: This routine uses LAPACK routines xTREXC.
544: Level: developer
546: .seealso: EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
547: @*/
548: PetscErrorCode EPSSortDenseSchur(EPS eps,PetscInt n_,PetscInt k,PetscScalar *T,PetscInt ldt_,PetscScalar *Q,PetscScalar *wr,PetscScalar *wi)
549: {
550: #if defined(SLEPC_MISSING_LAPACK_TREXC)
552: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
553: #else
555: PetscScalar re,im;
556: PetscInt i,j,pos,result;
557: PetscBLASInt ifst,ilst,info,n,ldt;
558: #if !defined(PETSC_USE_COMPLEX)
559: PetscScalar *work;
560: #endif
561:
563: PetscLogEventBegin(EPS_Dense,0,0,0,0);
564: n = PetscBLASIntCast(n_);
565: ldt = PetscBLASIntCast(ldt_);
566: #if !defined(PETSC_USE_COMPLEX)
567: PetscMalloc(n*sizeof(PetscScalar),&work);
568: #endif
569:
570: /* selection sort */
571: for (i=k;i<n-1;i++) {
572: re = wr[i];
573: im = wi[i];
574: pos = 0;
575: j=i+1; /* j points to the next eigenvalue */
576: #if !defined(PETSC_USE_COMPLEX)
577: if (im != 0) j=i+2;
578: #endif
579: /* find minimum eigenvalue */
580: for (;j<n;j++) {
581: EPSCompareEigenvalues(eps,re,im,wr[j],wi[j],&result);
582: if (result > 0) {
583: re = wr[j];
584: im = wi[j];
585: pos = j;
586: }
587: #if !defined(PETSC_USE_COMPLEX)
588: if (wi[j] != 0) j++;
589: #endif
590: }
591: if (pos) {
592: /* interchange blocks */
593: ifst = PetscBLASIntCast(pos + 1);
594: ilst = PetscBLASIntCast(i + 1);
595: #if !defined(PETSC_USE_COMPLEX)
596: LAPACKtrexc_("V",&n,T,&ldt,Q,&n,&ifst,&ilst,work,&info);
597: #else
598: LAPACKtrexc_("V",&n,T,&ldt,Q,&n,&ifst,&ilst,&info);
599: #endif
600: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
601: /* recover original eigenvalues from T matrix */
602: for (j=i;j<n;j++) {
603: wr[j] = T[j*ldt+j];
604: #if !defined(PETSC_USE_COMPLEX)
605: if (j<n-1 && T[j*ldt+j+1] != 0.0) {
606: /* complex conjugate eigenvalue */
607: wi[j] = PetscSqrtReal(PetscAbsReal(T[j*ldt+j+1])) *
608: PetscSqrtReal(PetscAbsReal(T[(j+1)*ldt+j]));
609: wr[j+1] = wr[j];
610: wi[j+1] = -wi[j];
611: j++;
612: } else
613: #endif
614: wi[j] = 0.0;
615: }
616: }
617: #if !defined(PETSC_USE_COMPLEX)
618: if (wi[i] != 0) i++;
619: #endif
620: }
622: #if !defined(PETSC_USE_COMPLEX)
623: PetscFree(work);
624: #endif
625: PetscLogEventEnd(EPS_Dense,0,0,0,0);
626: return(0);
628: #endif
629: }
633: /*@
634: EPSSortDenseSchurGeneralized - Reorders a generalized Schur decomposition.
636: Not Collective
638: Input Parameters:
639: + eps - the eigensolver context
640: . n - dimension of the matrix
641: . k0 - first active column
642: . k1 - last column to be ordered
643: - ldt - leading dimension of T
645: Input/Output Parameters:
646: + T,S - the upper (quasi-)triangular matrices
647: . Q,Z - the orthogonal matrix of Schur vectors
648: . wr - pointer to the array to store the computed eigenvalues
649: - wi - imaginary part of the eigenvalues (only when using real numbers)
651: Notes:
652: This function reorders the eigenvalues in wr,wi located in positions k0
653: to n according to the sort order specified in EPSetWhicheigenpairs.
654: The selection sort is the method used to sort the eigenvalues, and it
655: stops when the column k1-1 is ordered. The Schur decomposition Z*T*Z^T,
656: is also reordered by means of rotations so that eigenvalues in the
657: diagonal blocks of T follow the same order.
659: T,S,Q and Z are overwritten.
660:
661: This routine uses LAPACK routines xTGEXC.
663: Level: developer
665: .seealso: EPSSortDenseSchur(), EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()
666: @*/
667: PetscErrorCode EPSSortDenseSchurGeneralized(EPS eps,PetscInt n_,PetscInt k0,PetscInt k1,PetscScalar *T,PetscScalar *S,PetscInt ldt_,PetscScalar *Q,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
668: {
669: #if defined(SLEPC_MISSING_LAPACK_TGEXC) || !defined(PETSC_USE_COMPLEX) && (defined(SLEPC_MISSING_LAPACK_LAMCH) || defined(SLEPC_MISSING_LAPACK_LAG2))
671: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TGEXC/LAMCH/LAG2 - Lapack routines are unavailable.");
672: #else
674: PetscScalar re,im;
675: PetscInt i,j,result,pos;
676: PetscBLASInt ione = 1,ifst,ilst,info,n,ldt;
677: #if !defined(PETSC_USE_COMPLEX)
678: PetscBLASInt lwork;
679: PetscScalar *work,safmin,scale1,scale2,tmp;
680: #endif
681:
683: PetscLogEventBegin(EPS_Dense,0,0,0,0);
684: n = PetscBLASIntCast(n_);
685: ldt = PetscBLASIntCast(ldt_);
686: #if !defined(PETSC_USE_COMPLEX)
687: lwork = -1;
688: LAPACKtgexc_(&ione,&ione,&n,T,&ldt,S,&ldt,Q,&n,Z,&n,&ione,&ione,&tmp,&lwork,&info);
689: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTGEXC %d",info);
690: lwork = (PetscBLASInt)tmp;
691: safmin = LAPACKlamch_("S");
692: PetscMalloc(lwork*sizeof(PetscScalar),&work);
693: #endif
694:
695: /* selection sort */
696: for (i=k0;i<PetscMin(n-1,k1);i++) {
697: re = wr[i];
698: im = wi[i];
699: pos = 0;
700: j=i+1; /* j points to the next eigenvalue */
701: #if !defined(PETSC_USE_COMPLEX)
702: if (im != 0) j=i+2;
703: #endif
704: /* find minimum eigenvalue */
705: for (;j<n;j++) {
706: EPSCompareEigenvalues(eps,re,im,wr[j],wi[j],&result);
707: if (result > 0) {
708: re = wr[j];
709: im = wi[j];
710: pos = j;
711: }
712: #if !defined(PETSC_USE_COMPLEX)
713: if (wi[j] != 0) j++;
714: #endif
715: }
716: if (pos) {
717: /* interchange blocks */
718: ifst = PetscBLASIntCast(pos + 1);
719: ilst = PetscBLASIntCast(i + 1);
720: #if !defined(PETSC_USE_COMPLEX)
721: LAPACKtgexc_(&ione,&ione,&n,T,&ldt,S,&ldt,Q,&n,Z,&n,&ifst,&ilst,work,&lwork,&info);
722: #else
723: LAPACKtgexc_(&ione,&ione,&n,T,&ldt,S,&ldt,Q,&n,Z,&n,&ifst,&ilst,&info);
724: #endif
725: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTGEXC %d",info);
726: /* recover original eigenvalues from T and S matrices */
727: for (j=k0;j<n;j++) {
728: #if !defined(PETSC_USE_COMPLEX)
729: if (j<n-1 && T[j*ldt+j+1] != 0.0) {
730: /* complex conjugate eigenvalue */
731: LAPACKlag2_(T+j*ldt+j,&ldt,S+j*ldt+j,&ldt,&safmin,&scale1,&scale2,&re,&tmp,&im);
732: wr[j] = re / scale1;
733: wi[j] = im / scale1;
734: wr[j+1] = tmp / scale2;
735: wi[j+1] = -wi[j];
736: j++;
737: } else
738: #endif
739: {
740: if (S[j*ldt+j] == 0.0) {
741: if (PetscRealPart(T[j*ldt+j]) < 0.0) wr[j] = PETSC_MIN_REAL;
742: else wr[j] = PETSC_MAX_REAL;
743: } else wr[j] = T[j*ldt+j] / S[j*ldt+j];
744: wi[j] = 0.0;
745: }
746: }
747: }
748: #if !defined(PETSC_USE_COMPLEX)
749: if (wi[i] != 0) i++;
750: #endif
751: }
753: #if !defined(PETSC_USE_COMPLEX)
754: PetscFree(work);
755: #endif
756: PetscLogEventEnd(EPS_Dense,0,0,0,0);
757: return(0);
759: #endif
760: }
764: /*@
765: EPSDenseTridiagonal - Solves a real tridiagonal Hermitian Eigenvalue Problem.
767: Not Collective
769: Input Parameters:
770: + n - dimension of the eigenproblem
771: . D - pointer to the array containing the diagonal elements
772: - E - pointer to the array containing the off-diagonal elements
774: Output Parameters:
775: + w - pointer to the array to store the computed eigenvalues
776: - V - pointer to the array to store the eigenvectors
778: Notes:
779: If V is PETSC_NULL then the eigenvectors are not computed.
781: This routine use LAPACK routines xSTEVR.
783: Level: developer
785: .seealso: EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
786: @*/
787: PetscErrorCode EPSDenseTridiagonal(PetscInt n_,PetscReal *D,PetscReal *E,PetscReal *w,PetscScalar *V)
788: {
789: #if defined(SLEPC_MISSING_LAPACK_STEVR)
791: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"STEVR - Lapack routine is unavailable.");
792: #else
794: PetscReal abstol = 0.0,vl,vu,*work;
795: PetscBLASInt il,iu,m,*isuppz,n,lwork,*iwork,liwork,info;
796: const char *jobz;
797: #if defined(PETSC_USE_COMPLEX)
798: PetscInt i,j;
799: PetscReal *VV;
800: #endif
801:
803: PetscLogEventBegin(EPS_Dense,0,0,0,0);
804: n = PetscBLASIntCast(n_);
805: lwork = PetscBLASIntCast(20*n_);
806: liwork = PetscBLASIntCast(10*n_);
807: if (V) {
808: jobz = "V";
809: #if defined(PETSC_USE_COMPLEX)
810: PetscMalloc(n*n*sizeof(PetscReal),&VV);
811: #endif
812: } else jobz = "N";
813: PetscMalloc(2*n*sizeof(PetscBLASInt),&isuppz);
814: PetscMalloc(lwork*sizeof(PetscReal),&work);
815: PetscMalloc(liwork*sizeof(PetscBLASInt),&iwork);
816: #if defined(PETSC_USE_COMPLEX)
817: LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,VV,&n,isuppz,work,&lwork,iwork,&liwork,&info);
818: #else
819: LAPACKstevr_(jobz,"A",&n,D,E,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info);
820: #endif
821: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack DSTEVR %d",info);
822: #if defined(PETSC_USE_COMPLEX)
823: if (V) {
824: for (i=0;i<n;i++)
825: for (j=0;j<n;j++)
826: V[i*n+j] = VV[i*n+j];
827: PetscFree(VV);
828: }
829: #endif
830: PetscFree(isuppz);
831: PetscFree(work);
832: PetscFree(iwork);
833: PetscLogEventEnd(EPS_Dense,0,0,0,0);
834: return(0);
835: #endif
836: }
840: /*
841: DenseSelectedEvec - Computes a selected eigenvector of matrix in Schur form.
843: Input Parameters:
844: S - (quasi-)triangular matrix (dimension nv, leading dimension lds)
845: U - orthogonal transformation matrix (dimension nv, leading dimension nv)
846: i - which eigenvector to process
847: iscomplex - true if a complex conjugate pair (in real scalars)
849: Output parameters:
850: Y - computed eigenvector, 2 columns if iscomplex=true (leading dimension nv)
852: Workspace:
853: work is workspace to store 3*nv scalars, nv booleans and nv reals
854: */
855: PetscErrorCode DenseSelectedEvec(PetscScalar *S,PetscInt lds_,PetscScalar *U,PetscScalar *Y,PetscInt i,PetscBool iscomplex,PetscInt nv_,PetscScalar *work)
856: {
857: #if defined(SLEPC_MISSING_LAPACK_TREVC)
859: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable.");
860: #else
862: PetscInt k;
863: PetscBLASInt mm,mout,info,lds,nv,inc = 1;
864: PetscScalar tmp,done=1.0,zero=0.0;
865: PetscReal norm;
866: PetscBool *select=(PetscBool*)(work+4*nv_);
867: #if defined(PETSC_USE_COMPLEX)
868: PetscReal *rwork=(PetscReal*)(work+3*nv_);
869: #endif
872: lds = PetscBLASIntCast(lds_);
873: nv = PetscBLASIntCast(nv_);
874: for (k=0;k<nv;k++) select[k] = PETSC_FALSE;
876: /* Compute eigenvectors Y of S */
877: mm = iscomplex? 2: 1;
878: select[i] = PETSC_TRUE;
879: #if !defined(PETSC_USE_COMPLEX)
880: if (iscomplex) select[i+1] = PETSC_TRUE;
881: LAPACKtrevc_("R","S",select,&nv,S,&lds,PETSC_NULL,&nv,Y,&nv,&mm,&mout,work,&info);
882: #else
883: LAPACKtrevc_("R","S",select,&nv,S,&lds,PETSC_NULL,&nv,Y,&nv,&mm,&mout,work,rwork,&info);
884: #endif
885: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xTREVC %d",info);
886: if (mout != mm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Inconsistent arguments");
887: PetscMemcpy(work,Y,mout*nv*sizeof(PetscScalar));
889: /* accumulate and normalize eigenvectors */
890: BLASgemv_("N",&nv,&nv,&done,U,&nv,work,&inc,&zero,Y,&inc);
891: #if !defined(PETSC_USE_COMPLEX)
892: if (iscomplex) BLASgemv_("N",&nv,&nv,&done,U,&nv,work+nv,&inc,&zero,Y+nv,&inc);
893: #endif
894: mm = mm*nv;
895: norm = BLASnrm2_(&mm,Y,&inc);
896: tmp = 1.0 / norm;
897: BLASscal_(&mm,&tmp,Y,&inc);
898: return(0);
899: #endif
900: }