Actual source code: epsdefault.c
slepc-3.7.4 2017-05-17
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-2016, 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 && eps->purify) {
53: /* Purify eigenvectors */
54: BVCreateVec(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: if (eps->purify) {
96: BVCreateVec(eps->V,&v);
97: for (i=0;i<eps->nconv;i++) {
98: BVCopyVec(eps->V,i,v);
99: BVGetColumn(eps->V,i,&z);
100: STApply(eps->st,v,z);
101: BVRestoreColumn(eps->V,i,&z);
102: }
103: VecDestroy(&v);
104: }
106: /* normalization */
107: for (i=0;i<n;i++) {
108: #if !defined(PETSC_USE_COMPLEX)
109: if (eps->eigi[i] != 0.0) {
110: BVGetColumn(eps->V,i,&v);
111: BVGetColumn(eps->V,i+1,&v1);
112: VecNorm(v,NORM_2,&norm);
113: VecNorm(v1,NORM_2,&normi);
114: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
115: VecScale(v,tmp);
116: VecScale(v1,tmp);
117: BVRestoreColumn(eps->V,i,&v);
118: BVRestoreColumn(eps->V,i+1,&v1);
119: i++;
120: } else
121: #endif
122: {
123: BVGetColumn(eps->V,i,&v);
124: VecNormalize(v,NULL);
125: BVRestoreColumn(eps->V,i,&v);
126: }
127: }
128: return(0);
129: }
133: /*
134: EPSComputeVectors_Schur - Compute eigenvectors from the vectors
135: provided by the eigensolver. This version is intended for solvers
136: that provide Schur vectors. Given the partial Schur decomposition
137: OP*V=V*T, the following steps are performed:
138: 1) compute eigenvectors of T: T*Z=Z*D
139: 2) compute eigenvectors of OP: X=V*Z
140: */
141: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
142: {
144: PetscInt n,i;
145: Mat Z;
146: Vec w,z,v;
147: #if !defined(PETSC_USE_COMPLEX)
148: Vec v1;
149: PetscScalar tmp;
150: PetscReal norm,normi;
151: #endif
154: if (eps->ishermitian) {
155: if (eps->isgeneralized && !eps->ispositive) {
156: EPSComputeVectors_Indefinite(eps);
157: } else {
158: EPSComputeVectors_Hermitian(eps);
159: }
160: return(0);
161: }
162: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
164: /* right eigenvectors */
165: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
167: /* V = V * Z */
168: DSGetMat(eps->ds,DS_MAT_X,&Z);
169: BVSetActiveColumns(eps->V,0,n);
170: BVMultInPlace(eps->V,Z,0,n);
171: MatDestroy(&Z);
173: /* Purify eigenvectors */
174: if (eps->ispositive && eps->purify) {
175: BVCreateVec(eps->V,&w);
176: for (i=0;i<n;i++) {
177: BVCopyVec(eps->V,i,w);
178: BVGetColumn(eps->V,i,&z);
179: STApply(eps->st,w,z);
180: BVRestoreColumn(eps->V,i,&z);
181: }
182: VecDestroy(&w);
183: }
185: /* Fix eigenvectors if balancing was used */
186: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
187: for (i=0;i<n;i++) {
188: BVGetColumn(eps->V,i,&z);
189: VecPointwiseDivide(z,z,eps->D);
190: BVRestoreColumn(eps->V,i,&z);
191: }
192: }
194: /* normalize eigenvectors (when using purification or balancing) */
195: if ((eps->ispositive && eps->purify) || (eps->balance!=EPS_BALANCE_NONE && eps->D)) {
196: for (i=0;i<n;i++) {
197: #if !defined(PETSC_USE_COMPLEX)
198: if (eps->eigi[i] != 0.0) {
199: BVGetColumn(eps->V,i,&v);
200: BVGetColumn(eps->V,i+1,&v1);
201: VecNorm(v,NORM_2,&norm);
202: VecNorm(v1,NORM_2,&normi);
203: tmp = 1.0 / SlepcAbsEigenvalue(norm,normi);
204: VecScale(v,tmp);
205: VecScale(v1,tmp);
206: BVRestoreColumn(eps->V,i,&v);
207: BVRestoreColumn(eps->V,i+1,&v1);
208: i++;
209: } else
210: #endif
211: {
212: BVGetColumn(eps->V,i,&v);
213: VecNormalize(v,NULL);
214: BVRestoreColumn(eps->V,i,&v);
215: }
216: }
217: }
218: return(0);
219: }
223: /*@
224: EPSSetWorkVecs - Sets a number of work vectors into an EPS object.
226: Collective on EPS
228: Input Parameters:
229: + eps - eigensolver context
230: - nw - number of work vectors to allocate
232: Developers Note:
233: This is PETSC_EXTERN because it may be required by user plugin EPS
234: implementations.
236: Level: developer
237: @*/
238: PetscErrorCode EPSSetWorkVecs(EPS eps,PetscInt nw)
239: {
241: Vec t;
244: if (eps->nwork < nw) {
245: VecDestroyVecs(eps->nwork,&eps->work);
246: eps->nwork = nw;
247: BVGetColumn(eps->V,0,&t);
248: VecDuplicateVecs(t,nw,&eps->work);
249: BVRestoreColumn(eps->V,0,&t);
250: PetscLogObjectParents(eps,nw,eps->work);
251: }
252: return(0);
253: }
257: /*
258: EPSSetWhichEigenpairs_Default - Sets the default value for which,
259: depending on the ST.
260: */
261: PetscErrorCode EPSSetWhichEigenpairs_Default(EPS eps)
262: {
263: PetscBool target;
267: PetscObjectTypeCompareAny((PetscObject)eps->st,&target,STSINVERT,STCAYLEY,"");
268: if (target) eps->which = EPS_TARGET_MAGNITUDE;
269: else eps->which = EPS_LARGEST_MAGNITUDE;
270: return(0);
271: }
275: /*
276: EPSConvergedRelative - Checks convergence relative to the eigenvalue.
277: */
278: PetscErrorCode EPSConvergedRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
279: {
280: PetscReal w;
283: w = SlepcAbsEigenvalue(eigr,eigi);
284: *errest = res/w;
285: return(0);
286: }
290: /*
291: EPSConvergedAbsolute - Checks convergence absolutely.
292: */
293: PetscErrorCode EPSConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
294: {
296: *errest = res;
297: return(0);
298: }
302: /*
303: EPSConvergedNorm - Checks convergence relative to the eigenvalue and
304: the matrix norms.
305: */
306: PetscErrorCode EPSConvergedNorm(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
307: {
308: PetscReal w;
311: w = SlepcAbsEigenvalue(eigr,eigi);
312: *errest = res / (eps->nrma + w*eps->nrmb);
313: return(0);
314: }
318: /*@C
319: EPSStoppingBasic - Default routine to determine whether the outer eigensolver
320: iteration must be stopped.
322: Collective on EPS
324: Input Parameters:
325: + eps - eigensolver context obtained from EPSCreate()
326: . its - current number of iterations
327: . max_it - maximum number of iterations
328: . nconv - number of currently converged eigenpairs
329: . nev - number of requested eigenpairs
330: - ctx - context (not used here)
332: Output Parameter:
333: . reason - result of the stopping test
335: Notes:
336: A positive value of reason indicates that the iteration has finished successfully
337: (converged), and a negative value indicates an error condition (diverged). If
338: the iteration needs to be continued, reason must be set to EPS_CONVERGED_ITERATING
339: (zero).
341: EPSStoppingBasic() will stop if all requested eigenvalues are converged, or if
342: the maximum number of iterations has been reached.
344: Use EPSSetStoppingTest() to provide your own test instead of using this one.
346: Level: advanced
348: .seealso: EPSSetStoppingTest(), EPSConvergedReason, EPSGetConvergedReason()
349: @*/
350: PetscErrorCode EPSStoppingBasic(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
351: {
355: *reason = EPS_CONVERGED_ITERATING;
356: if (nconv >= nev) {
357: PetscInfo2(eps,"Linear eigensolver finished successfully: %D eigenpairs converged at iteration %D\n",nconv,its);
358: *reason = EPS_CONVERGED_TOL;
359: } else if (its >= max_it) {
360: *reason = EPS_DIVERGED_ITS;
361: PetscInfo1(eps,"Linear eigensolver iteration reached maximum number of iterations (%D)\n",its);
362: }
363: return(0);
364: }
368: /*
369: EPSComputeRitzVector - Computes the current Ritz vector.
371: Simple case (complex scalars or real scalars with Zi=NULL):
372: x = V*Zr (V is a basis of nv vectors, Zr has length nv)
374: Split case:
375: x = V*Zr y = V*Zi (Zr and Zi have length nv)
376: */
377: PetscErrorCode EPSComputeRitzVector(EPS eps,PetscScalar *Zr,PetscScalar *Zi,BV V,Vec x,Vec y)
378: {
380: PetscReal norm;
381: #if !defined(PETSC_USE_COMPLEX)
382: Vec z;
383: #endif
386: /* compute eigenvector */
387: BVMultVec(V,1.0,0.0,x,Zr);
389: /* purify eigenvector in positive generalized problems */
390: if (eps->ispositive && eps->purify) {
391: STApply(eps->st,x,y);
392: if (eps->ishermitian) {
393: BVNormVec(eps->V,y,NORM_2,&norm);
394: } else {
395: VecNorm(y,NORM_2,&norm);
396: }
397: VecScale(y,1.0/norm);
398: VecCopy(y,x);
399: }
400: /* fix eigenvector if balancing is used */
401: if (!eps->ishermitian && eps->balance!=EPS_BALANCE_NONE && eps->D) {
402: VecPointwiseDivide(x,x,eps->D);
403: VecNormalize(x,&norm);
404: }
405: #if !defined(PETSC_USE_COMPLEX)
406: /* compute imaginary part of eigenvector */
407: if (Zi) {
408: BVMultVec(V,1.0,0.0,y,Zi);
409: if (eps->ispositive) {
410: BVCreateVec(V,&z);
411: STApply(eps->st,y,z);
412: VecNorm(z,NORM_2,&norm);
413: VecScale(z,1.0/norm);
414: VecCopy(z,y);
415: VecDestroy(&z);
416: }
417: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
418: VecPointwiseDivide(y,y,eps->D);
419: VecNormalize(y,&norm);
420: }
421: } else
422: #endif
423: { VecSet(y,0.0); }
424: return(0);
425: }
429: /*
430: EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
431: diagonal matrix to be applied for balancing in non-Hermitian problems.
432: */
433: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
434: {
435: Vec z,p,r;
436: PetscInt i,j;
437: PetscReal norma;
438: PetscScalar *pz,*pD;
439: const PetscScalar *pr,*pp;
440: PetscRandom rand;
441: PetscErrorCode ierr;
444: EPSSetWorkVecs(eps,3);
445: BVGetRandomContext(eps->V,&rand);
446: r = eps->work[0];
447: p = eps->work[1];
448: z = eps->work[2];
449: VecSet(eps->D,1.0);
451: for (j=0;j<eps->balance_its;j++) {
453: /* Build a random vector of +-1's */
454: VecSetRandom(z,rand);
455: VecGetArray(z,&pz);
456: for (i=0;i<eps->nloc;i++) {
457: if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
458: else pz[i]=1.0;
459: }
460: VecRestoreArray(z,&pz);
462: /* Compute p=DA(D\z) */
463: VecPointwiseDivide(r,z,eps->D);
464: STApply(eps->st,r,p);
465: VecPointwiseMult(p,p,eps->D);
466: if (j==0) {
467: /* Estimate the matrix inf-norm */
468: VecAbs(p);
469: VecMax(p,NULL,&norma);
470: }
471: if (eps->balance == EPS_BALANCE_TWOSIDE) {
472: /* Compute r=D\(A'Dz) */
473: VecPointwiseMult(z,z,eps->D);
474: STApplyTranspose(eps->st,z,r);
475: VecPointwiseDivide(r,r,eps->D);
476: }
478: /* Adjust values of D */
479: VecGetArrayRead(r,&pr);
480: VecGetArrayRead(p,&pp);
481: VecGetArray(eps->D,&pD);
482: for (i=0;i<eps->nloc;i++) {
483: if (eps->balance == EPS_BALANCE_TWOSIDE) {
484: if (PetscAbsScalar(pp[i])>eps->balance_cutoff*norma && pr[i]!=0.0)
485: pD[i] *= PetscSqrtReal(PetscAbsScalar(pr[i]/pp[i]));
486: } else {
487: if (pp[i]!=0.0) pD[i] *= 1.0/PetscAbsScalar(pp[i]);
488: }
489: }
490: VecRestoreArrayRead(r,&pr);
491: VecRestoreArrayRead(p,&pp);
492: VecRestoreArray(eps->D,&pD);
493: }
494: return(0);
495: }