Actual source code: default.c
1: /*
2: This file contains some simple default routines for common operations.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
25: #include <slepcblaslapack.h>
29: PetscErrorCode EPSReset_Default(EPS eps)
30: {
34: EPSDefaultFreeWork(eps);
35: EPSFreeSolution(eps);
36: return(0);
37: }
41: PetscErrorCode EPSBackTransform_Default(EPS eps)
42: {
46: STBackTransform(eps->OP,eps->nconv,eps->eigr,eps->eigi);
47: return(0);
48: }
52: /*
53: EPSComputeVectors_Default - Compute eigenvectors from the vectors
54: provided by the eigensolver. This version just copies the vectors
55: and is intended for solvers such as power that provide the eigenvector.
56: */
57: PetscErrorCode EPSComputeVectors_Default(EPS eps)
58: {
60: eps->evecsavailable = PETSC_TRUE;
61: return(0);
62: }
66: /*
67: EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
68: using purification for generalized eigenproblems.
69: */
70: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
71: {
73: PetscInt i;
74: PetscReal norm;
75: Vec w;
78: if (eps->isgeneralized) {
79: /* Purify eigenvectors */
80: VecDuplicate(eps->V[0],&w);
81: for (i=0;i<eps->nconv;i++) {
82: VecCopy(eps->V[i],w);
83: STApply(eps->OP,w,eps->V[i]);
84: IPNorm(eps->ip,eps->V[i],&norm);
85: VecScale(eps->V[i],1.0/norm);
86: }
87: VecDestroy(&w);
88: }
89: eps->evecsavailable = PETSC_TRUE;
90: return(0);
91: }
95: /*
96: EPSComputeVectors_Indefinite - similar to the Schur version but
97: for indefinite problems
98: */
99: PetscErrorCode EPSComputeVectors_Indefinite(EPS eps)
100: {
102: PetscInt n,ld,i;
103: PetscScalar *Z;
104: Vec v;
105: #if !defined(PETSC_USE_COMPLEX)
106: PetscScalar tmp;
107: PetscReal norm,normi;
108: #endif
111: DSGetLeadingDimension(eps->ds,&ld);
112: DSGetDimensions(eps->ds,&n,PETSC_NULL,PETSC_NULL,PETSC_NULL);
113: DSVectors(eps->ds,DS_MAT_X,PETSC_NULL,PETSC_NULL);
114: DSGetArray(eps->ds,DS_MAT_X,&Z);
115: SlepcUpdateVectors(n,eps->V,0,n,Z,ld,PETSC_FALSE);
116: DSRestoreArray(eps->ds,DS_MAT_X,&Z);
117: /* purification */
118: VecDuplicate(eps->V[0],&v);
119: for (i=0;i<eps->nconv;i++) {
120: VecCopy(eps->V[i],v);
121: STApply(eps->OP,v,eps->V[i]);
122: }
123: VecDestroy(&v);
124: /* normalization */
125: for (i=0;i<n;i++) {
126: #if !defined(PETSC_USE_COMPLEX)
127: if (eps->eigi[i] != 0.0) {
128: VecNorm(eps->V[i],NORM_2,&norm);
129: VecNorm(eps->V[i+1],NORM_2,&normi);
130: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
131: VecScale(eps->V[i],tmp);
132: VecScale(eps->V[i+1],tmp);
133: i++;
134: } else
135: #endif
136: {
137: VecNormalize(eps->V[i],PETSC_NULL);
138: }
139: }
140: eps->evecsavailable = PETSC_TRUE;
141: return(0);
142: }
146: /*
147: EPSComputeVectors_Schur - Compute eigenvectors from the vectors
148: provided by the eigensolver. This version is intended for solvers
149: that provide Schur vectors. Given the partial Schur decomposition
150: OP*V=V*T, the following steps are performed:
151: 1) compute eigenvectors of T: T*Z=Z*D
152: 2) compute eigenvectors of OP: X=V*Z
153: If left eigenvectors are required then also do Z'*T=D*Z', Y=W*Z
154: */
155: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
156: {
158: PetscInt n,i,ld;
159: PetscScalar *Z;
160: #if !defined(PETSC_USE_COMPLEX)
161: PetscScalar tmp;
162: PetscReal norm,normi;
163: #endif
164: Vec w;
165:
167: if (eps->ishermitian) {
168: if(eps->isgeneralized && !eps->ispositive) {
169: EPSComputeVectors_Indefinite(eps);
170: } else {
171: EPSComputeVectors_Hermitian(eps);
172: }
173: return(0);
174: }
175: DSGetLeadingDimension(eps->ds,&ld);
176: DSGetDimensions(eps->ds,&n,PETSC_NULL,PETSC_NULL,PETSC_NULL);
178: /* right eigenvectors */
179: DSVectors(eps->ds,DS_MAT_X,PETSC_NULL,PETSC_NULL);
181: /* V = V * Z */
182: DSGetArray(eps->ds,DS_MAT_X,&Z);
183: SlepcUpdateVectors(n,eps->V,0,n,Z,ld,PETSC_FALSE);
184: DSRestoreArray(eps->ds,DS_MAT_X,&Z);
186: /* Purify eigenvectors */
187: if (eps->ispositive) {
188: VecDuplicate(eps->V[0],&w);
189: for (i=0;i<n;i++) {
190: VecCopy(eps->V[i],w);
191: STApply(eps->OP,w,eps->V[i]);
192: }
193: VecDestroy(&w);
194: }
196: /* Fix eigenvectors if balancing was used */
197: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
198: for (i=0;i<n;i++) {
199: VecPointwiseDivide(eps->V[i],eps->V[i],eps->D);
200: }
201: }
203: /* normalize eigenvectors (when using purification or balancing) */
204: if (eps->ispositive || (eps->balance!=EPS_BALANCE_NONE && eps->D)) {
205: for (i=0;i<n;i++) {
206: #if !defined(PETSC_USE_COMPLEX)
207: if (eps->eigi[i] != 0.0) {
208: VecNorm(eps->V[i],NORM_2,&norm);
209: VecNorm(eps->V[i+1],NORM_2,&normi);
210: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
211: VecScale(eps->V[i],tmp);
212: VecScale(eps->V[i+1],tmp);
213: i++;
214: } else
215: #endif
216: {
217: VecNormalize(eps->V[i],PETSC_NULL);
218: }
219: }
220: }
221:
222: /* left eigenvectors */
223: if (eps->leftvecs) {
224: DSVectors(eps->ds,DS_MAT_Y,PETSC_NULL,PETSC_NULL);
225: /* W = W * Z */
226: DSGetArray(eps->ds,DS_MAT_Y,&Z);
227: SlepcUpdateVectors(n,eps->W,0,n,Z,ld,PETSC_FALSE);
228: DSRestoreArray(eps->ds,DS_MAT_Y,&Z);
229: }
230: eps->evecsavailable = PETSC_TRUE;
231: return(0);
232: }
236: /*
237: EPSDefaultGetWork - Gets a number of work vectors.
238: */
239: PetscErrorCode EPSDefaultGetWork(EPS eps,PetscInt nw)
240: {
244: if (eps->nwork != nw) {
245: VecDestroyVecs(eps->nwork,&eps->work);
246: eps->nwork = nw;
247: VecDuplicateVecs(eps->t,nw,&eps->work);
248: PetscLogObjectParents(eps,nw,eps->work);
249: }
250: return(0);
251: }
255: /*
256: EPSDefaultFreeWork - Free work vectors.
257: */
258: PetscErrorCode EPSDefaultFreeWork(EPS eps)
259: {
263: VecDestroyVecs(eps->nwork,&eps->work);
264: eps->nwork = 0;
265: return(0);
266: }
270: /*
271: EPSDefaultSetWhich - Sets the default value for which, depending on the ST.
272: */
273: PetscErrorCode EPSDefaultSetWhich(EPS eps)
274: {
275: PetscBool target;
279: PetscObjectTypeCompareAny((PetscObject)eps->OP,&target,STSINVERT,STCAYLEY,STFOLD,"");
280: if (target) eps->which = EPS_TARGET_MAGNITUDE;
281: else eps->which = EPS_LARGEST_MAGNITUDE;
282: return(0);
283: }
287: /*
288: EPSConvergedEigRelative - Checks convergence relative to the eigenvalue.
289: */
290: PetscErrorCode EPSConvergedEigRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
291: {
292: PetscReal w;
295: w = SlepcAbsEigenvalue(eigr,eigi);
296: *errest = res/w;
297: return(0);
298: }
302: /*
303: EPSConvergedAbsolute - Checks convergence absolutely.
304: */
305: PetscErrorCode EPSConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
306: {
308: *errest = res;
309: return(0);
310: }
314: /*
315: EPSConvergedNormRelative - Checks convergence relative to the eigenvalue and
316: the matrix norms.
317: */
318: PetscErrorCode EPSConvergedNormRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
319: {
320: PetscReal w;
323: w = SlepcAbsEigenvalue(eigr,eigi);
324: *errest = res / (eps->nrma + w*eps->nrmb);
325: return(0);
326: }
330: /*
331: EPSComputeRitzVector - Computes the current Ritz vector.
333: Simple case (complex scalars or real scalars with Zi=NULL):
334: x = V*Zr (V is an array of nv vectors, Zr has length nv)
336: Split case:
337: x = V*Zr y = V*Zi (Zr and Zi have length nv)
338: */
339: PetscErrorCode EPSComputeRitzVector(EPS eps,PetscScalar *Zr,PetscScalar *Zi,Vec *V,PetscInt nv,Vec x,Vec y)
340: {
342: PetscReal norm;
343: #if !defined(PETSC_USE_COMPLEX)
344: Vec z;
345: #endif
346:
348: /* compute eigenvector */
349: SlepcVecMAXPBY(x,0.0,1.0,nv,Zr,V);
351: /* purify eigenvector in positive generalized problems */
352: if (eps->ispositive) {
353: STApply(eps->OP,x,y);
354: if (eps->ishermitian) {
355: IPNorm(eps->ip,y,&norm);
356: } else {
357: VecNorm(y,NORM_2,&norm);
358: }
359: VecScale(y,1.0/norm);
360: VecCopy(y,x);
361: }
362: /* fix eigenvector if balancing is used */
363: if (!eps->ishermitian && eps->balance!=EPS_BALANCE_NONE && eps->D) {
364: VecPointwiseDivide(x,x,eps->D);
365: VecNormalize(x,&norm);
366: }
367: #if !defined(PETSC_USE_COMPLEX)
368: /* compute imaginary part of eigenvector */
369: if (Zi) {
370: SlepcVecMAXPBY(y,0.0,1.0,nv,Zi,V);
371: if (eps->ispositive) {
372: VecDuplicate(V[0],&z);
373: STApply(eps->OP,y,z);
374: VecNorm(z,NORM_2,&norm);
375: VecScale(z,1.0/norm);
376: VecCopy(z,y);
377: VecDestroy(&z);
378: }
379: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
380: VecPointwiseDivide(y,y,eps->D);
381: VecNormalize(y,&norm);
382: }
383: } else
384: #endif
385: { VecSet(y,0.0); }
386: return(0);
387: }
391: /*
392: EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
393: diagonal matrix to be applied for balancing in non-Hermitian problems.
394: */
395: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
396: {
397: Vec z,p,r;
398: PetscInt i,j;
399: PetscReal norma;
400: PetscScalar *pz,*pD;
401: const PetscScalar *pr,*pp;
402: PetscErrorCode ierr;
405: VecDuplicate(eps->V[0],&r);
406: VecDuplicate(eps->V[0],&p);
407: VecDuplicate(eps->V[0],&z);
408: VecSet(eps->D,1.0);
410: for (j=0;j<eps->balance_its;j++) {
412: /* Build a random vector of +-1's */
413: SlepcVecSetRandom(z,eps->rand);
414: VecGetArray(z,&pz);
415: for (i=0;i<eps->nloc;i++) {
416: if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
417: else pz[i]=1.0;
418: }
419: VecRestoreArray(z,&pz);
421: /* Compute p=DA(D\z) */
422: VecPointwiseDivide(r,z,eps->D);
423: STApply(eps->OP,r,p);
424: VecPointwiseMult(p,p,eps->D);
425: if (j==0) {
426: /* Estimate the matrix inf-norm */
427: VecAbs(p);
428: VecMax(p,PETSC_NULL,&norma);
429: }
430: if (eps->balance == EPS_BALANCE_TWOSIDE) {
431: /* Compute r=D\(A'Dz) */
432: VecPointwiseMult(z,z,eps->D);
433: STApplyTranspose(eps->OP,z,r);
434: VecPointwiseDivide(r,r,eps->D);
435: }
436:
437: /* Adjust values of D */
438: VecGetArrayRead(r,&pr);
439: VecGetArrayRead(p,&pp);
440: VecGetArray(eps->D,&pD);
441: for (i=0;i<eps->nloc;i++) {
442: if (eps->balance == EPS_BALANCE_TWOSIDE) {
443: if (PetscAbsScalar(pp[i])>eps->balance_cutoff*norma && pr[i]!=0.0)
444: pD[i] *= PetscSqrtReal(PetscAbsScalar(pr[i]/pp[i]));
445: } else {
446: if (pp[i]!=0.0) pD[i] *= 1.0/PetscAbsScalar(pp[i]);
447: }
448: }
449: VecRestoreArrayRead(r,&pr);
450: VecRestoreArrayRead(p,&pp);
451: VecRestoreArray(eps->D,&pD);
452: }
454: VecDestroy(&r);
455: VecDestroy(&p);
456: VecDestroy(&z);
457: return(0);
458: }