Actual source code: dense.c
2: /*
3: This file implements interfaces to direct solvers in LAPACK
4: */
5: #include slepceps.h
6: #include slepcblaslapack.h
8: #define SWAP(a,b,t) {t=a;a=b;b=t;}
12: /*@
13: EPSSortEigenvalues - Sorts a list of eigenvalues according to a certain
14: criterion.
16: Not Collective
18: Input Parameters:
19: + n - number of eigenvalue in the list
20: . eig - pointer to the array containing the eigenvalues
21: . eigi - imaginary part of the eigenvalues (only when using real numbers)
22: . which - sorting criterion
23: - nev - number of wanted eigenvalues
25: Output Parameter:
26: . permout - resulting permutation
28: Notes:
29: The result is a list of indices in the original eigenvalue array
30: corresponding to the first nev eigenvalues sorted in the specified
31: criterion
33: Level: developer
35: .seealso: EPSDenseNHEPSorted(), EPSSetWhichEigenpairs()
36: @*/
37: int EPSSortEigenvalues(int n,PetscScalar *eig,PetscScalar *eigi,EPSWhich which,int nev,int *permout)
38: {
39: int ierr,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: #if defined(PETSC_USE_COMPLEX)
51: for (i=0; i<n; i++) { values[i] = PetscAbsScalar(eig[i]); }
52: #else
53: for (i=0; i<n; i++) { values[i] = LAlapy2_(&eig[i],&eigi[i]); }
54: #endif
55: break;
56: case EPS_LARGEST_REAL:
57: case EPS_SMALLEST_REAL:
58: for (i=0; i<n; i++) { values[i] = PetscRealPart(eig[i]); }
59: break;
60: case EPS_LARGEST_IMAGINARY:
61: case EPS_SMALLEST_IMAGINARY:
62: #if defined(PETSC_USE_COMPLEX)
63: for (i=0; i<n; i++) { values[i] = PetscImaginaryPart(eig[i]); }
64: #else
65: for (i=0; i<n; i++) { values[i] = PetscAbsReal(eigi[i]); }
66: #endif
67: break;
68: default: SETERRQ(1,"Wrong value of which");
69: }
71: PetscSortRealWithPermutation(n,values,perm);
73: switch(which) {
74: case EPS_LARGEST_MAGNITUDE:
75: case EPS_LARGEST_REAL:
76: case EPS_LARGEST_IMAGINARY:
77: for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
78: break;
79: case EPS_SMALLEST_MAGNITUDE:
80: case EPS_SMALLEST_REAL:
81: case EPS_SMALLEST_IMAGINARY:
82: for (i=0; i<nev; i++) { permout[i] = perm[i]; }
83: break;
84: default: SETERRQ(1,"Wrong value of which");
85: }
87: #if !defined(PETSC_USE_COMPLEX)
88: for (i=0; i<nev-1; i++) {
89: if (eigi[permout[i]] != 0.0) {
90: if (eig[permout[i]] == eig[permout[i+1]] &&
91: eigi[permout[i]] == -eigi[permout[i+1]] &&
92: eigi[permout[i]] < 0.0) {
93: int tmp;
94: SWAP(permout[i], permout[i+1], tmp);
95: }
96: i++;
97: }
98: }
99: #endif
101: PetscFree(values);
102: PetscFree(perm);
103: return(0);
104: }
108: /*@
109: EPSDenseNHEP - Solves a dense non-Hermitian Eigenvalue Problem.
111: Not Collective
113: Input Parameter:
114: + n - dimension of the eigenproblem
115: - A - pointer to the array containing the matrix values
117: Output Parameters:
118: + w - pointer to the array to store the computed eigenvalues
119: . wi - imaginary part of the eigenvalues (only when using real numbers)
120: - V - pointer to the array to store the eigenvectors
122: Notes:
123: If V is PETSC_NULL then the eigenvectors are not computed.
125: Matrix A is overwritten.
126:
127: This routine uses LAPACK routines xGEEV.
129: Level: developer
131: .seealso: EPSDenseNHEPSorted()
132: @*/
133: int EPSDenseNHEP(int n,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V)
134: {
136:
137: #if defined(PETSC_HAVE_ESSL)
139: /* ESSL has a different calling sequence for dgeev() and zgeev() than
140: standard LAPACK */
141: PetscScalar *cwork;
142: PetscReal *work;
143: int i,clen,idummy,lwork,iopt;
146: #if !defined(PETSC_USE_COMPLEX)
147: clen = n;
148: #else
149: clen = 2*n;
150: #endif
151: PetscMalloc(clen*sizeof(PetscScalar),&cwork);
152: idummy = n;
153: lwork = 3*n;
154: PetscMalloc(lwork*sizeof(PetscReal),&work);
155: if (V) iopt = 1;
156: else iopt = 0;
157: LAgeev_(&iopt,A,&n,cwork,V,&n,&idummy,&n,work,&lwork);
158: PetscFree(work);
159: #if !defined(PETSC_USE_COMPLEX)
160: for (i=0; i<n; i++) {
161: w[i] = cwork[2*i];
162: wi[i] = cwork[2*i+1];
163: }
164: #else
165: for (i=0; i<n; i++) w[i] = cwork[i];
166: #endif
167: PetscFree(cwork);
169: #elif !defined(PETSC_USE_COMPLEX)
171: PetscScalar *work,sdummy;
172: int lwork;
173: char *jobvr;
176: lwork = 5*n;
177: PetscMalloc(lwork*sizeof(PetscScalar),&work);
178: if (V) jobvr = "V";
179: else jobvr = "N";
180: LAgeev_("N",jobvr,&n,A,&n,w,wi,&sdummy,&n,V,&n,work,&lwork,&ierr);
181: if (ierr) SETERRQ1(PETSC_ERR_LIB,"Error in LAPACK routine %d",ierr);
182: PetscFree(work);
184: #else
186: PetscScalar *work,sdummy;
187: PetscReal *rwork;
188: int lwork;
189: char *jobvr;
192: lwork = 5*n;
193: PetscMalloc(lwork*sizeof(PetscScalar),&work);
194: PetscMalloc(2*n*sizeof(PetscReal),&rwork);
195: if (V) jobvr = "V";
196: else jobvr = "N";
197: LAgeev_("N",jobvr,&n,A,&n,w,&sdummy,&n,V,&n,work,&lwork,rwork,&ierr);
198: if (ierr) SETERRQ1(PETSC_ERR_LIB,"Error in LAPACK routine %d",ierr);
199: PetscFree(work);
200: PetscFree(rwork);
202: #endif
204: return(0);
205: }
209: /*@
210: EPSDenseNHEPSorted - Solves a dense non-Hermitian Eigenvalue Problem and
211: then sorts the computed eigenpairs.
213: Not Collective
215: Input Parameter:
216: + n - dimension of the eigenproblem
217: - A - pointer to the array containing the matrix values
219: Output Parameters:
220: + w - pointer to the array to store the computed eigenvalues
221: . wi - imaginary part of the eigenvalues (only when using real numbers)
222: - V - pointer to the array to store the eigenvectors
224: Notes:
225: If V is PETSC_NULL then the eigenvectors are not computed.
227: Matrix A is overwritten.
229: Level: developer
231: .seealso: EPSDenseNHEP(), EPSSortEigenvalues()
232: @*/
233: int EPSDenseNHEPSorted(int n,PetscScalar *A,PetscScalar *w,PetscScalar *wi,PetscScalar *V,int m,EPSWhich which)
234: {
235: int i,ierr,*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: }