Actual source code: epsdefault.c
slepc-3.5.2 2014-10-10
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-2014, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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 <slepcvec.h>
29: PetscErrorCode EPSBackTransform_Default(EPS eps)
30: {
34: STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi);
35: return(0);
36: }
40: /*
41: EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
42: using purification for generalized eigenproblems.
43: */
44: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
45: {
47: PetscInt i;
48: PetscReal norm;
49: Vec w,z;
52: if (eps->isgeneralized) {
53: /* Purify eigenvectors */
54: BVGetVec(eps->V,&w);
55: for (i=0;i<eps->nconv;i++) {
56: BVCopyVec(eps->V,i,w);
57: BVGetColumn(eps->V,i,&z);
58: STApply(eps->st,w,z);
59: BVRestoreColumn(eps->V,i,&z);
60: BVNormColumn(eps->V,i,NORM_2,&norm);
61: BVScaleColumn(eps->V,i,1.0/norm);
62: }
63: VecDestroy(&w);
64: }
65: return(0);
66: }
70: /*
71: EPSComputeVectors_Indefinite - similar to the Schur version but
72: for indefinite problems
73: */
74: PetscErrorCode EPSComputeVectors_Indefinite(EPS eps)
75: {
77: PetscInt n,i;
78: Mat X;
79: Vec v,z;
80: #if !defined(PETSC_USE_COMPLEX)
81: Vec v1;
82: PetscScalar tmp;
83: PetscReal norm,normi;
84: #endif
87: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
88: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
89: DSGetMat(eps->ds,DS_MAT_X,&X);
90: BVSetActiveColumns(eps->V,0,n);
91: BVMultInPlace(eps->V,X,0,n);
92: MatDestroy(&X);
94: /* purification */
95: BVGetVec(eps->V,&v);
96: for (i=0;i<eps->nconv;i++) {
97: BVCopyVec(eps->V,i,v);
98: BVGetColumn(eps->V,i,&z);
99: STApply(eps->st,v,z);
100: BVRestoreColumn(eps->V,i,&z);
101: }
102: VecDestroy(&v);
104: /* normalization */
105: for (i=0;i<n;i++) {
106: #if !defined(PETSC_USE_COMPLEX)
107: if (eps->eigi[i] != 0.0) {
108: BVGetColumn(eps->V,i,&v);
109: BVGetColumn(eps->V,i+1,&v1);
110: VecNorm(v,NORM_2,&norm);
111: VecNorm(v1,NORM_2,&normi);
112: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
113: VecScale(v,tmp);
114: VecScale(v1,tmp);
115: BVRestoreColumn(eps->V,i,&v);
116: BVRestoreColumn(eps->V,i+1,&v1);
117: i++;
118: } else
119: #endif
120: {
121: BVGetColumn(eps->V,i,&v);
122: VecNormalize(v,NULL);
123: BVRestoreColumn(eps->V,i,&v);
124: }
125: }
126: return(0);
127: }
131: /*
132: EPSComputeVectors_Schur - Compute eigenvectors from the vectors
133: provided by the eigensolver. This version is intended for solvers
134: that provide Schur vectors. Given the partial Schur decomposition
135: OP*V=V*T, the following steps are performed:
136: 1) compute eigenvectors of T: T*Z=Z*D
137: 2) compute eigenvectors of OP: X=V*Z
138: */
139: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
140: {
142: PetscInt n,i;
143: Mat Z;
144: Vec w,z,v;
145: #if !defined(PETSC_USE_COMPLEX)
146: Vec v1;
147: PetscScalar tmp;
148: PetscReal norm,normi;
149: #endif
152: if (eps->ishermitian) {
153: if (eps->isgeneralized && !eps->ispositive) {
154: EPSComputeVectors_Indefinite(eps);
155: } else {
156: EPSComputeVectors_Hermitian(eps);
157: }
158: return(0);
159: }
160: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
162: /* right eigenvectors */
163: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
165: /* V = V * Z */
166: DSGetMat(eps->ds,DS_MAT_X,&Z);
167: BVSetActiveColumns(eps->V,0,n);
168: BVMultInPlace(eps->V,Z,0,n);
169: MatDestroy(&Z);
171: /* Purify eigenvectors */
172: if (eps->ispositive) {
173: BVGetVec(eps->V,&w);
174: for (i=0;i<n;i++) {
175: BVCopyVec(eps->V,i,w);
176: BVGetColumn(eps->V,i,&z);
177: STApply(eps->st,w,z);
178: BVRestoreColumn(eps->V,i,&z);
179: }
180: VecDestroy(&w);
181: }
183: /* Fix eigenvectors if balancing was used */
184: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
185: for (i=0;i<n;i++) {
186: BVGetColumn(eps->V,i,&z);
187: VecPointwiseDivide(z,z,eps->D);
188: BVRestoreColumn(eps->V,i,&z);
189: }
190: }
192: /* normalize eigenvectors (when using purification or balancing) */
193: if (eps->ispositive || (eps->balance!=EPS_BALANCE_NONE && eps->D)) {
194: for (i=0;i<n;i++) {
195: #if !defined(PETSC_USE_COMPLEX)
196: if (eps->eigi[i] != 0.0) {
197: BVGetColumn(eps->V,i,&v);
198: BVGetColumn(eps->V,i+1,&v1);
199: VecNorm(v,NORM_2,&norm);
200: VecNorm(v1,NORM_2,&normi);
201: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
202: VecScale(v,tmp);
203: VecScale(v1,tmp);
204: BVRestoreColumn(eps->V,i,&v);
205: BVRestoreColumn(eps->V,i+1,&v1);
206: i++;
207: } else
208: #endif
209: {
210: BVGetColumn(eps->V,i,&v);
211: VecNormalize(v,NULL);
212: BVRestoreColumn(eps->V,i,&v);
213: }
214: }
215: }
216: return(0);
217: }
221: /*@
222: EPSSetWorkVecs - Sets a number of work vectors into a EPS object.
224: Collective on EPS
226: Input Parameters:
227: + eps - eigensolver context
228: - nw - number of work vectors to allocate
230: Developers Note:
231: This is PETSC_EXTERN because it may be required by user plugin EPS
232: implementations.
234: Level: developer
235: @*/
236: PetscErrorCode EPSSetWorkVecs(EPS eps,PetscInt nw)
237: {
239: Vec t;
242: if (eps->nwork != nw) {
243: VecDestroyVecs(eps->nwork,&eps->work);
244: eps->nwork = nw;
245: BVGetColumn(eps->V,0,&t);
246: VecDuplicateVecs(t,nw,&eps->work);
247: BVRestoreColumn(eps->V,0,&t);
248: PetscLogObjectParents(eps,nw,eps->work);
249: }
250: return(0);
251: }
255: /*
256: EPSSetWhichEigenpairs_Default - Sets the default value for which,
257: depending on the ST.
258: */
259: PetscErrorCode EPSSetWhichEigenpairs_Default(EPS eps)
260: {
261: PetscBool target;
265: PetscObjectTypeCompareAny((PetscObject)eps->st,&target,STSINVERT,STCAYLEY,"");
266: if (target) eps->which = EPS_TARGET_MAGNITUDE;
267: else eps->which = EPS_LARGEST_MAGNITUDE;
268: return(0);
269: }
273: /*
274: EPSConvergedEigRelative - Checks convergence relative to the eigenvalue.
275: */
276: PetscErrorCode EPSConvergedEigRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
277: {
278: PetscReal w;
281: w = SlepcAbsEigenvalue(eigr,eigi);
282: *errest = res/w;
283: return(0);
284: }
288: /*
289: EPSConvergedAbsolute - Checks convergence absolutely.
290: */
291: PetscErrorCode EPSConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
292: {
294: *errest = res;
295: return(0);
296: }
300: /*
301: EPSConvergedNormRelative - Checks convergence relative to the eigenvalue and
302: the matrix norms.
303: */
304: PetscErrorCode EPSConvergedNormRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
305: {
306: PetscReal w;
309: w = SlepcAbsEigenvalue(eigr,eigi);
310: *errest = res / (eps->nrma + w*eps->nrmb);
311: return(0);
312: }
316: /*
317: EPSComputeRitzVector - Computes the current Ritz vector.
319: Simple case (complex scalars or real scalars with Zi=NULL):
320: x = V*Zr (V is a basis of nv vectors, Zr has length nv)
322: Split case:
323: x = V*Zr y = V*Zi (Zr and Zi have length nv)
324: */
325: PetscErrorCode EPSComputeRitzVector(EPS eps,PetscScalar *Zr,PetscScalar *Zi,BV V,Vec x,Vec y)
326: {
328: PetscReal norm;
329: #if !defined(PETSC_USE_COMPLEX)
330: Vec z;
331: #endif
334: /* compute eigenvector */
335: BVMultVec(V,1.0,0.0,x,Zr);
337: /* purify eigenvector in positive generalized problems */
338: if (eps->ispositive) {
339: STApply(eps->st,x,y);
340: if (eps->ishermitian) {
341: BVNormVec(eps->V,y,NORM_2,&norm);
342: } else {
343: VecNorm(y,NORM_2,&norm);
344: }
345: VecScale(y,1.0/norm);
346: VecCopy(y,x);
347: }
348: /* fix eigenvector if balancing is used */
349: if (!eps->ishermitian && eps->balance!=EPS_BALANCE_NONE && eps->D) {
350: VecPointwiseDivide(x,x,eps->D);
351: VecNormalize(x,&norm);
352: }
353: #if !defined(PETSC_USE_COMPLEX)
354: /* compute imaginary part of eigenvector */
355: if (Zi) {
356: BVMultVec(V,1.0,0.0,y,Zi);
357: if (eps->ispositive) {
358: BVGetVec(V,&z);
359: STApply(eps->st,y,z);
360: VecNorm(z,NORM_2,&norm);
361: VecScale(z,1.0/norm);
362: VecCopy(z,y);
363: VecDestroy(&z);
364: }
365: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
366: VecPointwiseDivide(y,y,eps->D);
367: VecNormalize(y,&norm);
368: }
369: } else
370: #endif
371: { VecSet(y,0.0); }
372: return(0);
373: }
377: /*
378: EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
379: diagonal matrix to be applied for balancing in non-Hermitian problems.
380: */
381: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
382: {
383: Vec z,p,r;
384: PetscInt i,j;
385: PetscReal norma;
386: PetscScalar *pz,*pD;
387: const PetscScalar *pr,*pp;
388: PetscErrorCode ierr;
391: BVGetVec(eps->V,&r);
392: BVGetVec(eps->V,&p);
393: BVGetVec(eps->V,&z);
394: VecSet(eps->D,1.0);
396: for (j=0;j<eps->balance_its;j++) {
398: /* Build a random vector of +-1's */
399: VecSetRandom(z,eps->rand);
400: VecGetArray(z,&pz);
401: for (i=0;i<eps->nloc;i++) {
402: if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
403: else pz[i]=1.0;
404: }
405: VecRestoreArray(z,&pz);
407: /* Compute p=DA(D\z) */
408: VecPointwiseDivide(r,z,eps->D);
409: STApply(eps->st,r,p);
410: VecPointwiseMult(p,p,eps->D);
411: if (j==0) {
412: /* Estimate the matrix inf-norm */
413: VecAbs(p);
414: VecMax(p,NULL,&norma);
415: }
416: if (eps->balance == EPS_BALANCE_TWOSIDE) {
417: /* Compute r=D\(A'Dz) */
418: VecPointwiseMult(z,z,eps->D);
419: STApplyTranspose(eps->st,z,r);
420: VecPointwiseDivide(r,r,eps->D);
421: }
423: /* Adjust values of D */
424: VecGetArrayRead(r,&pr);
425: VecGetArrayRead(p,&pp);
426: VecGetArray(eps->D,&pD);
427: for (i=0;i<eps->nloc;i++) {
428: if (eps->balance == EPS_BALANCE_TWOSIDE) {
429: if (PetscAbsScalar(pp[i])>eps->balance_cutoff*norma && pr[i]!=0.0)
430: pD[i] *= PetscSqrtReal(PetscAbsScalar(pr[i]/pp[i]));
431: } else {
432: if (pp[i]!=0.0) pD[i] *= 1.0/PetscAbsScalar(pp[i]);
433: }
434: }
435: VecRestoreArrayRead(r,&pr);
436: VecRestoreArrayRead(p,&pp);
437: VecRestoreArray(eps->D,&pD);
438: }
440: VecDestroy(&r);
441: VecDestroy(&p);
442: VecDestroy(&z);
443: return(0);
444: }