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
175: Output Parameters:
176: + w - pointer to the array to store the computed eigenvalues
177: - V - pointer to the array to store the eigenvectors
179: Notes:
180: If V is PETSC_NULL then the eigenvectors are not computed.
182: Matrix A is overwritten.
183:
184: This routine uses LAPACK routines DSYEVR or ZHEEVR.
186: Level: developer
188: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseGHEP()
189: @*/
190: PetscErrorCode EPSDenseHEP(int n,PetscScalar *A,PetscReal *w,PetscScalar *V)
191: {
192: #if defined(SLEPC_MISSING_LAPACK_SYEVR) || defined(SLEPC_MISSING_LAPACK_HEEVR)
194: SETERRQ(PETSC_ERR_SUP,"DSYEVR/ZHEEVR - Lapack routine is unavailable.");
195: #else
197: PetscReal abstol = 0.0,vl,vu;
198: PetscScalar *work;
199: int il,iu,m,*isuppz,*iwork,liwork = 10*n,info;
200: const char *jobz;
201: #if defined(PETSC_USE_COMPLEX)
202: PetscReal *rwork;
203: int lwork = 18*n,lrwork = 24*n;
204: #else
205: int lwork = 26*n;
206: #endif
209: if (V) jobz = "V";
210: else jobz = "N";
211: PetscMalloc(2*n*sizeof(int),&isuppz);
212: PetscMalloc(lwork*sizeof(PetscScalar),&work);
213: PetscMalloc(liwork*sizeof(int),&iwork);
214: #if defined(PETSC_USE_COMPLEX)
215: PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
216: LAPACKsyevr_(jobz,"A","U",&n,A,&n,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,rwork,&lrwork,iwork,&liwork,&info,1,1,1);
217: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEEVR %d",info);
218: PetscFree(rwork);
219: #else
220: LAPACKsyevr_(jobz,"A","U",&n,A,&n,&vl,&vu,&il,&iu,&abstol,&m,w,V,&n,isuppz,work,&lwork,iwork,&liwork,&info,1,1,1);
221: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYEVR %d",info);
222: #endif
223: PetscFree(isuppz);
224: PetscFree(work);
225: PetscFree(iwork);
226: return(0);
227: #endif
228: }
232: /*@
233: EPSDenseGHEP - Solves a dense Generalized Hermitian Eigenvalue Problem.
235: Not Collective
237: Input Parameters:
238: + n - dimension of the eigenproblem
239: . A - pointer to the array containing the matrix values for A
240: - B - pointer to the array containing the matrix values for B
242: Output Parameters:
243: + w - pointer to the array to store the computed eigenvalues
244: - V - pointer to the array to store the eigenvectors
246: Notes:
247: If V is PETSC_NULL then the eigenvectors are not computed.
249: Matrices A and B are overwritten.
250:
251: This routine uses LAPACK routines DSYGVD or ZHEGVD.
253: Level: developer
255: .seealso: EPSDenseNHEP(), EPSDenseGNHEP(), EPSDenseHEP()
256: @*/
257: PetscErrorCode EPSDenseGHEP(int n,PetscScalar *A,PetscScalar *B,PetscReal *w,PetscScalar *V)
258: {
259: #if defined(SLEPC_MISSING_LAPACK_SYGVD) || defined(SLEPC_MISSING_LAPACK_HEGVD)
261: SETERRQ(PETSC_ERR_SUP,"DSYGVD/ZHEGVD - Lapack routine is unavailable.");
262: #else
264: PetscScalar *work;
265: int itype = 1,*iwork,info,
266: liwork = V ? 5*n+3 : 1;
267: const char *jobz;
268: #if defined(PETSC_USE_COMPLEX)
269: PetscReal *rwork;
270: int lwork = V ? n*n+2*n : n+1,
271: lrwork = V ? 2*n*n+5*n+1 : n;
272: #else
273: int lwork = V ? 2*n*n+6*n+1 : 2*n+1;
274: #endif
277: if (V) jobz = "V";
278: else jobz = "N";
279: PetscMalloc(lwork*sizeof(PetscScalar),&work);
280: PetscMalloc(liwork*sizeof(int),&iwork);
281: #if defined(PETSC_USE_COMPLEX)
282: PetscMalloc(lrwork*sizeof(PetscReal),&rwork);
283: LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,rwork,&lrwork,iwork,&liwork,&info,1,1);
284: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack ZHEGVD %d",info);
285: PetscFree(rwork);
286: #else
287: LAPACKsygvd_(&itype,jobz,"U",&n,A,&n,B,&n,w,work,&lwork,iwork,&liwork,&info,1,1);
288: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack DSYGVD %d",info);
289: #endif
290: if (V) {
291: PetscMemcpy(V,A,n*n*sizeof(PetscScalar));
292: }
293: PetscFree(work);
294: PetscFree(iwork);
295: return(0);
296: #endif
297: }
301: /*@
302: EPSDenseSchur - Computes the upper (quasi-)triangular form of a dense
303: upper Hessenberg matrix.
305: Not Collective
307: Input Parameters:
308: + n - dimension of the matrix
309: . k - first active column
310: - ldh - leading dimension of H
312: Input/Output Parameters:
313: + H - on entry, the upper Hessenber matrix; on exit, the upper
314: (quasi-)triangular matrix (T)
315: - Z - on entry, initial transformation matrix; on exit, orthogonal
316: matrix of Schur vectors
318: Output Parameters:
319: + wr - pointer to the array to store the computed eigenvalues
320: - wi - imaginary part of the eigenvalues (only when using real numbers)
322: Notes:
323: This function computes the (real) Schur decomposition of an upper
324: Hessenberg matrix H: H*Z = Z*T, where T is an upper (quasi-)triangular
325: matrix (returned in H), and Z is the orthogonal matrix of Schur vectors.
326: Eigenvalues are extracted from the diagonal blocks of T and returned in
327: wr,wi. Transformations are accumulated in Z so that on entry it can
328: contain the transformation matrix associated to the Hessenberg reduction.
330: Only active columns (from k to n) are computed.
332: Both H and Z are overwritten.
333:
334: This routine uses LAPACK routines xHSEQR.
336: Level: developer
338: .seealso: EPSSortDenseSchur()
339: @*/
340: PetscErrorCode EPSDenseSchur(int n,int k,PetscScalar *H,int ldh,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
341: {
342: #if defined(SLEPC_MISSING_LAPACK_HSEQR)
344: SETERRQ(PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable.");
345: #else
347: int ilo,lwork,info;
348: PetscScalar *work;
349: #if !defined(PETSC_USE_COMPLEX)
350: int j;
351: #endif
352:
354: lwork = n;
355: PetscMalloc(lwork*sizeof(PetscScalar),&work);
356: ilo = k+1;
357: #if !defined(PETSC_USE_COMPLEX)
358: LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,wi,Z,&n,work,&lwork,&info,1,1);
359: for (j=0;j<k;j++) {
360: if (j==n-1 || H[j*ldh+j+1] == 0.0) {
361: /* real eigenvalue */
362: wr[j] = H[j*ldh+j];
363: wi[j] = 0.0;
364: } else {
365: /* complex eigenvalue */
366: wr[j] = H[j*ldh+j];
367: wr[j+1] = H[j*ldh+j];
368: wi[j] = sqrt(PetscAbsReal(H[j*ldh+j+1])) *
369: sqrt(PetscAbsReal(H[(j+1)*ldh+j]));
370: wi[j+1] = -wi[j];
371: j++;
372: }
373: }
374: #else
375: LAPACKhseqr_("S","V",&n,&ilo,&n,H,&ldh,wr,Z,&n,work,&lwork,&info,1,1);
376: #endif
377: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);
379: PetscFree(work);
380: return(0);
381: #endif
382: }
386: /*@
387: EPSSortDenseSchur - Reorders the Schur decomposition computed by
388: EPSDenseSchur().
390: Not Collective
392: Input Parameters:
393: + n - dimension of the matrix
394: . k - first active column
395: - ldt - leading dimension of T
397: Input/Output Parameters:
398: + T - the upper (quasi-)triangular matrix
399: - Z - the orthogonal matrix of Schur vectors
401: Output Parameters:
402: + wr - pointer to the array to store the computed eigenvalues
403: - wi - imaginary part of the eigenvalues (only when using real numbers)
405: Notes:
406: This function reorders the eigenvalues in wr,wi located in positions k
407: to n in descending order of magnitude. The Schur decomposition Z*T*Z^T,
408: is also reordered by means of rotations so that eigenvalues in the
409: diagonal blocks of T follow the same order.
411: Both T and Z are overwritten.
412:
413: This routine uses LAPACK routines xTREXC.
415: Level: developer
417: .seealso: EPSDenseSchur()
418: @*/
419: PetscErrorCode EPSSortDenseSchur(int n,int k,PetscScalar *T,int ldt,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
420: {
421: #if defined(SLEPC_MISSING_LAPACK_TREXC)
423: SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
424: #else
425: int i,j,ifst,ilst,info,maxpos;
426: #if !defined(PETSC_USE_COMPLEX)
427: PetscScalar *work;
429: #endif
430: PetscReal max,m;
431:
434: #if !defined(PETSC_USE_COMPLEX)
436: PetscMalloc(n*sizeof(PetscScalar),&work);
437:
438: for (i=k;i<n-1;i++) {
439: max = SlepcAbsEigenvalue(wr[i],wi[i]);
440: maxpos = 0;
441: for (j=i+1;j<n;j++) {
442: m = SlepcAbsEigenvalue(wr[j],wi[j]);
443: if (m > max) {
444: max = m;
445: maxpos = j;
446: }
447: if (wi[j] != 0) j++;
448: }
449: if (maxpos) {
450: ifst = maxpos + 1;
451: ilst = i + 1;
452: LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,work,&info,1);
453: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
454:
455: for (j=k;j<n;j++) {
456: if (j==n-1 || T[j*ldt+j+1] == 0.0) {
457: /* real eigenvalue */
458: wr[j] = T[j*ldt+j];
459: wi[j] = 0.0;
460: } else {
461: /* complex eigenvalue */
462: wr[j] = T[j*ldt+j];
463: wr[j+1] = T[j*ldt+j];
464: wi[j] = sqrt(PetscAbsReal(T[j*ldt+j+1])) *
465: sqrt(PetscAbsReal(T[(j+1)*ldt+j]));
466: wi[j+1] = -wi[j];
467: j++;
468: }
469: }
470: }
471: if (wi[i] != 0) i++;
472: }
473:
474: PetscFree(work);
476: #else /* PETSC_USE_COMPLEX */
478: for (i=k;i<n-1;i++) {
479: max = SlepcAbsEigenvalue(wr[i],wi[i]);
480: maxpos = 0;
481: for (j=i+1;j<n;j++) {
482: m = SlepcAbsEigenvalue(wr[j],wi[j]);
483: if (m > max) {
484: max = m;
485: maxpos = j;
486: }
487: }
488: if (maxpos) {
489: ifst = maxpos + 1;
490: ilst = i + 1;
491: LAPACKtrexc_("V",&n,T,&ldt,Z,&n,&ifst,&ilst,&info,1);
492: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
494: for (j=k;j<n;j++) {
495: wr[j] = T[j*(ldt+1)];
496: }
497: }
498: }
500: #endif
501:
502: return(0);
504: #endif
505: }