Actual source code: dense.c
1: /*
2: This file contains routines for handling small size dense problems.
3: */
4: #include slepceps.h
5: #include slepcblaslapack.h
7: #define SWAP(a,b,t) {t=a;a=b;b=t;}
11: /*@
12: EPSSortEigenvalues - Sorts a list of eigenvalues according to a certain
13: criterion.
15: Not Collective
17: Input Parameters:
18: + n - number of eigenvalue in the list
19: . eig - pointer to the array containing the eigenvalues
20: . eigi - imaginary part of the eigenvalues (only when using real numbers)
21: . which - sorting criterion
22: - nev - number of wanted eigenvalues
24: Output Parameter:
25: . permout - resulting permutation
27: Notes:
28: The result is a list of indices in the original eigenvalue array
29: corresponding to the first nev eigenvalues sorted in the specified
30: criterion
32: Level: developer
34: .seealso: EPSDenseNHEPSorted(), EPSSetWhichEigenpairs()
35: @*/
36: PetscErrorCode EPSSortEigenvalues(int n,PetscScalar *eig,PetscScalar *eigi,EPSWhich which,int nev,int *permout)
37: {
38: PetscErrorCode ierr;
39: int i,*perm;
40: PetscReal *values;
43: PetscMalloc(n*sizeof(int),&perm);
44: PetscMalloc(n*sizeof(PetscReal),&values);
45: for (i=0; i<n; i++) { perm[i] = i;}
47: switch(which) {
48: case EPS_LARGEST_MAGNITUDE:
49: case EPS_SMALLEST_MAGNITUDE:
50: for (i=0; i<n; i++) { values[i] = SlepcAbsEigenvalue(eig[i],eigi[i]); }
51: break;
52: case EPS_LARGEST_REAL:
53: case EPS_SMALLEST_REAL:
54: for (i=0; i<n; i++) { values[i] = PetscRealPart(eig[i]); }
55: break;
56: case EPS_LARGEST_IMAGINARY:
57: case EPS_SMALLEST_IMAGINARY:
58: #if defined(PETSC_USE_COMPLEX)
59: for (i=0; i<n; i++) { values[i] = PetscImaginaryPart(eig[i]); }
60: #else
61: for (i=0; i<n; i++) { values[i] = PetscAbsReal(eigi[i]); }
62: #endif
63: break;
64: default: SETERRQ(1,"Wrong value of which");
65: }
67: PetscSortRealWithPermutation(n,values,perm);
69: switch(which) {
70: case EPS_LARGEST_MAGNITUDE:
71: case EPS_LARGEST_REAL:
72: case EPS_LARGEST_IMAGINARY:
73: for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
74: break;
75: case EPS_SMALLEST_MAGNITUDE:
76: case EPS_SMALLEST_REAL:
77: case EPS_SMALLEST_IMAGINARY:
78: for (i=0; i<nev; i++) { permout[i] = perm[i]; }
79: break;
80: default: SETERRQ(1,"Wrong value of which");
81: }
83: #if !defined(PETSC_USE_COMPLEX)
84: for (i=0; i<nev-1; i++) {
85: if (eigi[permout[i]] != 0.0) {
86: if (eig[permout[i]] == eig[permout[i+1]] &&
87: eigi[permout[i]] == -eigi[permout[i+1]] &&
88: eigi[permout[i]] < 0.0) {
89: int tmp;
90: SWAP(permout[i], permout[i+1], tmp);
91: }
92: i++;
93: }
94: }
95: #endif
97: PetscFree(values);
98: PetscFree(perm);
99: return(0);
100: }
104: /*@
105: EPSDenseNHEP - Solves a dense 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
113: Output Parameters:
114: + w - pointer to the array to store the computed eigenvalues
115: . wi - imaginary part of the eigenvalues (only when using real numbers)
116: - V - pointer to the array to store the eigenvectors
118: Notes:
119: If V is PETSC_NULL then the eigenvectors are not computed.
121: Matrix A is overwritten.
122:
123: This routine uses LAPACK routines xGEEV.
125: Level: developer
127: .seealso: EPSDenseNHEPSorted()
128: @*/
129: PetscErrorCode EPSDenseNHEP(int n,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V)
130: {
131: PetscErrorCode ierr;
132:
133: #if defined(PETSC_HAVE_ESSL)
135: /* ESSL has a different calling sequence for dgeev() and zgeev() than
136: standard LAPACK */
137: PetscScalar *cwork;
138: PetscReal *work;
139: int i,clen,idummy,lwork,iopt;
142: #if !defined(PETSC_USE_COMPLEX)
143: clen = n;
144: #else
145: clen = 2*n;
146: #endif
147: PetscMalloc(clen*sizeof(PetscScalar),&cwork);
148: idummy = n;
149: lwork = 3*n;
150: PetscMalloc(lwork*sizeof(PetscReal),&work);
151: if (V) iopt = 1;
152: else iopt = 0;
153: LAgeev_(&iopt,A,&n,cwork,V,&n,&idummy,&n,work,&lwork);
154: PetscFree(work);
155: #if !defined(PETSC_USE_COMPLEX)
156: for (i=0; i<n; i++) {
157: w[i] = cwork[2*i];
158: wi[i] = cwork[2*i+1];
159: }
160: #else
161: for (i=0; i<n; i++) w[i] = cwork[i];
162: #endif
163: PetscFree(cwork);
165: #elif !defined(PETSC_USE_COMPLEX)
167: PetscScalar *work,sdummy;
168: int lwork,info;
169: char *jobvr;
172: lwork = 5*n;
173: PetscMalloc(lwork*sizeof(PetscScalar),&work);
174: if (V) jobvr = "V";
175: else jobvr = "N";
176: LAgeev_("N",jobvr,&n,A,&n,w,wi,&sdummy,&n,V,&n,work,&lwork,&info);
177: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGEEV %d",info);
178: PetscFree(work);
180: #else
182: PetscScalar *work,sdummy;
183: PetscReal *rwork;
184: int lwork,info;
185: char *jobvr;
188: #if defined(PETSC_MISSING_LAPACK_GEEV)
189: SETERRQ(PETSC_ERR_SUP,"GEEV - Lapack routine is unavailable.");
190: #endif
191: lwork = 5*n;
192: PetscMalloc(lwork*sizeof(PetscScalar),&work);
193: PetscMalloc(2*n*sizeof(PetscReal),&rwork);
194: if (V) jobvr = "V";
195: else jobvr = "N";
196: LAgeev_("N",jobvr,&n,A,&n,w,&sdummy,&n,V,&n,work,&lwork,rwork,&info);
197: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGEEV %d",info);
198: PetscFree(work);
199: PetscFree(rwork);
201: #endif
203: return(0);
204: }
208: /*@
209: EPSDenseNHEPSorted - Solves a dense non-Hermitian Eigenvalue Problem and
210: then sorts the computed eigenpairs.
212: Not Collective
214: Input Parameters:
215: + n - dimension of the eigenproblem
216: - A - pointer to the array containing the matrix values
218: Output Parameters:
219: + w - pointer to the array to store the computed eigenvalues
220: . wi - imaginary part of the eigenvalues (only when using real numbers)
221: - V - pointer to the array to store the eigenvectors
223: Notes:
224: If V is PETSC_NULL then the eigenvectors are not computed.
226: Matrix A is overwritten.
228: Level: developer
230: .seealso: EPSDenseNHEP(), EPSSortEigenvalues()
231: @*/
232: PetscErrorCode EPSDenseNHEPSorted(int n,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V,int m,EPSWhich which)
233: {
234: PetscErrorCode ierr;
235: int i,*perm,iwork[100];
236: PetscScalar *realpart,*imagpart,*vectors,work[200];
239: if (m<=100) perm = iwork;
240: else { PetscMalloc(m*sizeof(int),&perm); }
241: if (n<=100) { realpart = work; imagpart = work+100; }
242: else {
243: PetscMalloc(n*sizeof(PetscScalar),&realpart);
244: PetscMalloc(n*sizeof(PetscScalar),&imagpart);
245: }
246: if (V) {
247: PetscMalloc(n*n*sizeof(PetscScalar),&vectors);
248: } else vectors = PETSC_NULL;
250: EPSDenseNHEP(n,A,realpart,imagpart,vectors);
252: EPSSortEigenvalues(n,realpart,imagpart,which,m,perm);
253: for (i=0; i<m; i++) {
254: w[i] = realpart[perm[i]];
255: #if !defined(PETSC_USE_COMPLEX)
256: wi[i] = imagpart[perm[i]];
257: #endif
258: if (V) {
259: PetscMemcpy(V+i*n,vectors+perm[i]*n,n*sizeof(PetscScalar));
260: }
261: }
263: if (m>100) { PetscFree(perm); }
264: if (n>100) {
265: PetscFree(realpart);
266: PetscFree(imagpart);
267: }
268: if (V) {
269: PetscFree(vectors);
270: }
272: return(0);
273: }
277: /*@
278: EPSDenseSchur - Computes the upper (quasi-)triangular form of a dense
279: upper Hessenberg matrix.
281: Not Collective
283: Input Parameters:
284: + n - dimension of the matrix
285: - k - first active column
287: Input/Output Parameters:
288: + H - on entry, the upper Hessenber matrix; on exit, the upper
289: (quasi-)triangular matrix (T)
290: - Z - on entry, initial transformation matrix; on exit, orthogonal
291: matrix of Schur vectors
293: Output Parameters:
294: + wr - pointer to the array to store the computed eigenvalues
295: - wi - imaginary part of the eigenvalues (only when using real numbers)
297: Notes:
298: This function computes the (real) Schur decomposition of an upper
299: Hessenberg matrix H: H*Z = Z*T, where T is an upper (quasi-)triangular
300: matrix (returned in H), and Z is the orthogonal matrix of Schur vectors.
301: Eigenvalues are extracted from the diagonal blocks of T and returned in
302: wr,wi. Transformations are accumulated in Z so that on entry it can
303: contain the transformation matrix associated to the Hessenberg reduction.
305: Only active columns (from k to n) are computed.
307: Both H and Z are overwritten.
308:
309: This routine uses LAPACK routines xHSEQR.
311: Level: developer
313: .seealso: EPSSortDenseSchur()
314: @*/
315: PetscErrorCode EPSDenseSchur(int n,int k,PetscScalar *H,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
316: {
317: PetscErrorCode ierr;
318: int ilo,lwork,info;
319: PetscScalar *work;
320: #if !defined(PETSC_USE_COMPLEX)
321: int j;
322: #endif
323:
325: #if defined(PETSC_BLASLAPACK_ESSL_ONLY)
326: SETERRQ(PETSC_ERR_SUP,"HSEQR - Lapack routine is unavailable.");
327: #endif
329: lwork = n;
330: PetscMalloc(lwork*sizeof(PetscScalar),&work);
331: ilo = k+1;
332: #if !defined(PETSC_USE_COMPLEX)
333: LAhseqr_("S","V",&n,&ilo,&n,H,&n,wr,wi,Z,&n,work,&lwork,&info,1,1);
334: for (j=0;j<k;j++) {
335: if (j==n-1 || H[j*n+j+1] == 0.0) {
336: /* real eigenvalue */
337: wr[j] = H[j*n+j];
338: wi[j] = 0.0;
339: } else {
340: /* complex eigenvalue */
341: wr[j] = H[j*n+j];
342: wr[j+1] = H[j*n+j];
343: wi[j] = sqrt(PetscAbsReal(H[j*n+j+1])) *
344: sqrt(PetscAbsReal(H[(j+1)*n+j]));
345: wi[j+1] = -wi[j];
346: j++;
347: }
348: }
349: #else
350: LAhseqr_("S","V",&n,&ilo,&n,H,&n,wr,Z,&n,work,&lwork,&info,1,1);
351: #endif
352: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xHSEQR %d",info);
354: PetscFree(work);
355: return(0);
356: }
360: /*@
361: EPSSortDenseSchur - Reorders the Schur decomposition computed by
362: EPSDenseSchur().
364: Not Collective
366: Input Parameters:
367: + n - dimension of the matrix
368: - k - first active column
370: Input/Output Parameters:
371: + T - the upper (quasi-)triangular matrix
372: - Z - the orthogonal matrix of Schur vectors
374: Output Parameters:
375: + wr - pointer to the array to store the computed eigenvalues
376: - wi - imaginary part of the eigenvalues (only when using real numbers)
378: Notes:
379: This function reorders the eigenvalues in wr,wi located in positions k
380: to n in ascending order of magnitude. The Schur decomposition Z*T*Z^T,
381: is also reordered by means of rotations so that eigenvalues in the
382: diagonal blocks of T follow the same order.
384: Both T and Z are overwritten.
385:
386: This routine uses LAPACK routines xTREXC.
388: Level: developer
390: .seealso: EPSDenseSchur()
391: @*/
392: PetscErrorCode EPSSortDenseSchur(int n,int k,PetscScalar *T,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi)
393: {
394: int i,j,ifst,ilst,info,maxpos;
395: #if !defined(PETSC_USE_COMPLEX)
396: PetscScalar *work;
397: PetscErrorCode ierr;
398: #endif
399: PetscReal max,m;
400:
402: #if defined(PETSC_BLASLAPACK_ESSL_ONLY)
403: SETERRQ(PETSC_ERR_SUP,"TREXC - Lapack routine is unavailable.");
404: #endif
405: #if !defined(PETSC_USE_COMPLEX)
406: PetscMalloc(n*sizeof(PetscScalar),&work);
407: #endif
409: for (i=k;i<n-1;i++) {
410: max = SlepcAbsEigenvalue(wr[i],wi[i]);
411: maxpos = 0;
412: for (j=i+1;j<n;j++) {
413: m = SlepcAbsEigenvalue(wr[j],wi[j]);
414: if (m > max) {
415: max = m;
416: maxpos = j;
417: }
418: #if !defined(PETSC_USE_COMPLEX)
419: if (wi[j] != 0) j++;
420: #endif
421: }
422: if (maxpos) {
423: ifst = maxpos + 1;
424: ilst = i + 1;
425: #if !defined(PETSC_USE_COMPLEX)
426: LAtrexc_("V",&n,T,&n,Z,&n,&ifst,&ilst,work,&info,1);
427: #else
428: LAtrexc_("V",&n,T,&n,Z,&n,&ifst,&ilst,&info,1);
429: #endif
430: if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xTREXC %d",info);
432: #if !defined(PETSC_USE_COMPLEX)
433: for (j=k;j<n;j++) {
434: if (j==n-1 || T[j*n+j+1] == 0.0) {
435: /* real eigenvalue */
436: wr[j] = T[j*n+j];
437: wi[j] = 0.0;
438: } else {
439: /* complex eigenvalue */
440: wr[j] = T[j*n+j];
441: wr[j+1] = T[j*n+j];
442: wi[j] = sqrt(PetscAbsReal(T[j*n+j+1])) *
443: sqrt(PetscAbsReal(T[(j+1)*n+j]));
444: wi[j+1] = -wi[j];
445: j++;
446: }
447: }
448: #else
449: for (j=k;j<n;j++) {
450: wr[j] = T[j*(n+1)];
451: }
452: #endif
453: }
454: #if !defined(PETSC_USE_COMPLEX)
455: if (wi[i] != 0) i++;
456: #endif
457: }
458:
459: #if !defined(PETSC_USE_COMPLEX)
460: PetscFree(work);
461: #endif
462: return(0);
463: }