Actual source code: epsdefault.c
slepc-main 2024-12-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: This file contains some simple default routines for common operations
12: */
14: #include <slepc/private/epsimpl.h>
15: #include <slepcvec.h>
17: PetscErrorCode EPSBackTransform_Default(EPS eps)
18: {
19: PetscFunctionBegin;
20: PetscCall(STBackTransform(eps->st,eps->nconv,eps->eigr,eps->eigi));
21: PetscFunctionReturn(PETSC_SUCCESS);
22: }
24: /*
25: EPSComputeVectors_Hermitian - Copies the Lanczos vectors as eigenvectors
26: using purification for generalized eigenproblems.
27: */
28: PetscErrorCode EPSComputeVectors_Hermitian(EPS eps)
29: {
30: PetscBool iscayley,indef;
31: Mat B,C;
33: PetscFunctionBegin;
34: if (eps->purify) {
35: PetscCall(EPS_Purify(eps,eps->nconv));
36: PetscCall(BVNormalize(eps->V,NULL));
37: } else {
38: /* In the case of Cayley transform, eigenvectors need to be B-normalized */
39: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley));
40: if (iscayley && eps->isgeneralized) {
41: PetscCall(STGetMatrix(eps->st,1,&B));
42: PetscCall(BVGetMatrix(eps->V,&C,&indef));
43: PetscCheck(!indef,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"The inner product should not be indefinite");
44: PetscCall(BVSetMatrix(eps->V,B,PETSC_FALSE));
45: PetscCall(BVNormalize(eps->V,NULL));
46: PetscCall(BVSetMatrix(eps->V,C,PETSC_FALSE)); /* restore original matrix */
47: }
48: }
49: PetscFunctionReturn(PETSC_SUCCESS);
50: }
52: /*
53: EPSComputeVectors_Indefinite - similar to the Schur version but
54: for indefinite problems
55: */
56: PetscErrorCode EPSComputeVectors_Indefinite(EPS eps)
57: {
58: PetscInt n;
59: Mat X;
61: PetscFunctionBegin;
62: PetscCall(DSGetDimensions(eps->ds,&n,NULL,NULL,NULL));
63: PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
64: PetscCall(DSGetMat(eps->ds,DS_MAT_X,&X));
65: PetscCall(BVMultInPlace(eps->V,X,0,n));
66: PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&X));
68: /* purification */
69: if (eps->purify) PetscCall(EPS_Purify(eps,eps->nconv));
71: /* normalization */
72: PetscCall(BVNormalize(eps->V,eps->eigi));
73: PetscFunctionReturn(PETSC_SUCCESS);
74: }
76: /*
77: EPSComputeVectors_Twosided - Adjust left eigenvectors in generalized problems: y = B^-* y.
78: */
79: PetscErrorCode EPSComputeVectors_Twosided(EPS eps)
80: {
81: PetscInt i;
82: Vec w,y;
84: PetscFunctionBegin;
85: if (!eps->twosided || !eps->isgeneralized) PetscFunctionReturn(PETSC_SUCCESS);
86: PetscCall(EPSSetWorkVecs(eps,1));
87: w = eps->work[0];
88: for (i=0;i<eps->nconv;i++) {
89: PetscCall(BVCopyVec(eps->W,i,w));
90: PetscCall(BVGetColumn(eps->W,i,&y));
91: PetscCall(STMatSolveHermitianTranspose(eps->st,w,y));
92: PetscCall(BVRestoreColumn(eps->W,i,&y));
93: }
94: PetscFunctionReturn(PETSC_SUCCESS);
95: }
97: /*
98: EPSComputeVectors_Schur - Compute eigenvectors from the vectors
99: provided by the eigensolver. This version is intended for solvers
100: that provide Schur vectors. Given the partial Schur decomposition
101: OP*V=V*T, the following steps are performed:
102: 1) compute eigenvectors of T: T*Z=Z*D
103: 2) compute eigenvectors of OP: X=V*Z
104: */
105: PetscErrorCode EPSComputeVectors_Schur(EPS eps)
106: {
107: PetscInt i;
108: Mat Z;
109: Vec z;
111: PetscFunctionBegin;
112: if (eps->ishermitian) {
113: if (eps->isgeneralized && !eps->ispositive) PetscCall(EPSComputeVectors_Indefinite(eps));
114: else PetscCall(EPSComputeVectors_Hermitian(eps));
115: PetscFunctionReturn(PETSC_SUCCESS);
116: }
118: /* right eigenvectors */
119: PetscCall(DSVectors(eps->ds,DS_MAT_X,NULL,NULL));
121: /* V = V * Z */
122: PetscCall(DSGetMat(eps->ds,DS_MAT_X,&Z));
123: PetscCall(BVMultInPlace(eps->V,Z,0,eps->nconv));
124: PetscCall(DSRestoreMat(eps->ds,DS_MAT_X,&Z));
126: /* Purify eigenvectors */
127: if (eps->purify) PetscCall(EPS_Purify(eps,eps->nconv));
129: /* Fix eigenvectors if balancing was used */
130: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
131: for (i=0;i<eps->nconv;i++) {
132: PetscCall(BVGetColumn(eps->V,i,&z));
133: PetscCall(VecPointwiseDivide(z,z,eps->D));
134: PetscCall(BVRestoreColumn(eps->V,i,&z));
135: }
136: }
138: /* normalize eigenvectors (when using purification or balancing) */
139: if (eps->purify || (eps->balance!=EPS_BALANCE_NONE && eps->D)) PetscCall(BVNormalize(eps->V,eps->eigi));
141: /* left eigenvectors */
142: if (eps->twosided) {
143: PetscCall(DSVectors(eps->ds,DS_MAT_Y,NULL,NULL));
144: /* W = W * Z */
145: PetscCall(DSGetMat(eps->ds,DS_MAT_Y,&Z));
146: PetscCall(BVMultInPlace(eps->W,Z,0,eps->nconv));
147: PetscCall(DSRestoreMat(eps->ds,DS_MAT_Y,&Z));
148: /* Fix left eigenvectors if balancing was used */
149: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
150: for (i=0;i<eps->nconv;i++) {
151: PetscCall(BVGetColumn(eps->W,i,&z));
152: PetscCall(VecPointwiseMult(z,z,eps->D));
153: PetscCall(BVRestoreColumn(eps->W,i,&z));
154: }
155: }
156: PetscCall(EPSComputeVectors_Twosided(eps));
157: /* normalize */
158: PetscCall(BVNormalize(eps->W,eps->eigi));
159: #if !defined(PETSC_USE_COMPLEX)
160: for (i=0;i<eps->nconv-1;i++) {
161: if (eps->eigi[i] != 0.0) {
162: if (eps->eigi[i] > 0.0) PetscCall(BVScaleColumn(eps->W,i+1,-1.0));
163: i++;
164: }
165: }
166: #endif
167: }
168: PetscFunctionReturn(PETSC_SUCCESS);
169: }
171: /*@
172: EPSSetWorkVecs - Sets a number of work vectors into an EPS object.
174: Collective
176: Input Parameters:
177: + eps - eigensolver context
178: - nw - number of work vectors to allocate
180: Developer Notes:
181: This is SLEPC_EXTERN because it may be required by user plugin EPS
182: implementations.
184: Level: developer
186: .seealso: EPSSetUp()
187: @*/
188: PetscErrorCode EPSSetWorkVecs(EPS eps,PetscInt nw)
189: {
190: Vec t;
192: PetscFunctionBegin;
195: PetscCheck(nw>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"nw must be > 0: nw = %" PetscInt_FMT,nw);
196: if (eps->nwork < nw) {
197: PetscCall(VecDestroyVecs(eps->nwork,&eps->work));
198: eps->nwork = nw;
199: PetscCall(BVGetColumn(eps->V,0,&t));
200: PetscCall(VecDuplicateVecs(t,nw,&eps->work));
201: PetscCall(BVRestoreColumn(eps->V,0,&t));
202: }
203: PetscFunctionReturn(PETSC_SUCCESS);
204: }
206: /*
207: EPSSetWhichEigenpairs_Default - Sets the default value for which,
208: depending on the ST.
209: */
210: PetscErrorCode EPSSetWhichEigenpairs_Default(EPS eps)
211: {
212: PetscBool target;
214: PetscFunctionBegin;
215: PetscCall(PetscObjectTypeCompareAny((PetscObject)eps->st,&target,STSINVERT,STCAYLEY,""));
216: if (target) eps->which = EPS_TARGET_MAGNITUDE;
217: else eps->which = EPS_LARGEST_MAGNITUDE;
218: PetscFunctionReturn(PETSC_SUCCESS);
219: }
221: /*
222: EPSConvergedRelative - Checks convergence relative to the eigenvalue.
223: */
224: PetscErrorCode EPSConvergedRelative(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
225: {
226: PetscReal w;
228: PetscFunctionBegin;
229: w = SlepcAbsEigenvalue(eigr,eigi);
230: *errest = (w!=0.0)? res/w: PETSC_MAX_REAL;
231: PetscFunctionReturn(PETSC_SUCCESS);
232: }
234: /*
235: EPSConvergedAbsolute - Checks convergence absolutely.
236: */
237: PetscErrorCode EPSConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
238: {
239: PetscFunctionBegin;
240: *errest = res;
241: PetscFunctionReturn(PETSC_SUCCESS);
242: }
244: /*
245: EPSConvergedNorm - Checks convergence relative to the eigenvalue and
246: the matrix norms.
247: */
248: PetscErrorCode EPSConvergedNorm(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
249: {
250: PetscReal w;
252: PetscFunctionBegin;
253: w = SlepcAbsEigenvalue(eigr,eigi);
254: *errest = res / (eps->nrma + w*eps->nrmb);
255: PetscFunctionReturn(PETSC_SUCCESS);
256: }
258: /*@C
259: EPSStoppingBasic - Default routine to determine whether the outer eigensolver
260: iteration must be stopped.
262: Collective
264: Input Parameters:
265: + eps - eigensolver context obtained from EPSCreate()
266: . its - current number of iterations
267: . max_it - maximum number of iterations
268: . nconv - number of currently converged eigenpairs
269: . nev - number of requested eigenpairs
270: - ctx - context (not used here)
272: Output Parameter:
273: . reason - result of the stopping test
275: Notes:
276: A positive value of reason indicates that the iteration has finished successfully
277: (converged), and a negative value indicates an error condition (diverged). If
278: the iteration needs to be continued, reason must be set to EPS_CONVERGED_ITERATING
279: (zero).
281: EPSStoppingBasic() will stop if all requested eigenvalues are converged, or if
282: the maximum number of iterations has been reached.
284: Use EPSSetStoppingTest() to provide your own test instead of using this one.
286: Level: advanced
288: .seealso: EPSSetStoppingTest(), EPSStoppingThreshold(), EPSConvergedReason, EPSGetConvergedReason()
289: @*/
290: PetscErrorCode EPSStoppingBasic(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
291: {
292: PetscFunctionBegin;
293: *reason = EPS_CONVERGED_ITERATING;
294: if (nconv >= nev) {
295: PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: %" PetscInt_FMT " eigenpairs converged at iteration %" PetscInt_FMT "\n",nconv,its));
296: *reason = EPS_CONVERGED_TOL;
297: } else if (its >= max_it) {
298: *reason = EPS_DIVERGED_ITS;
299: PetscCall(PetscInfo(eps,"Linear eigensolver iteration reached maximum number of iterations (%" PetscInt_FMT ")\n",its));
300: }
301: PetscFunctionReturn(PETSC_SUCCESS);
302: }
304: /*@C
305: EPSStoppingThreshold - Routine to determine whether the outer eigenvalue solver
306: iteration must be stopped, according to some threshold for the computed values.
308: Collective
310: Input Parameters:
311: + eps - eigenvalue solver context obtained from EPSCreate()
312: . its - current number of iterations
313: . max_it - maximum number of iterations
314: . nconv - number of currently converged eigenpairs (ignored here)
315: . nev - number of requested eigenpairs (ignored here)
316: - ctx - context containing additional data (EPSStoppingCtx)
318: Output Parameter:
319: . reason - result of the stopping test
321: Notes:
322: A positive value of reason indicates that the iteration has finished successfully
323: (converged), and a negative value indicates an error condition (diverged). If
324: the iteration needs to be continued, reason must be set to EPS_CONVERGED_ITERATING
325: (zero).
327: EPSStoppingThreshold() will stop when one of the computed eigenvalues is not
328: above/below the threshold given at EPSSetThreshold(). If a number of wanted
329: eigenvalues has been specified via EPSSetDimensions() then it is also taken into
330: account, and the solver will stop when one of the two conditions (threshold or
331: number of converged values) is met.
333: Use EPSSetStoppingTest() to provide your own test instead of using this one.
335: Level: advanced
337: .seealso: EPSSetStoppingTest(), EPSStoppingBasic(), EPSSetThreshold(), EPSSetDimensions(), EPSConvergedReason, EPSGetConvergedReason()
338: @*/
339: PetscErrorCode EPSStoppingThreshold(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
340: {
341: PetscReal thres,firstev,lastev;
342: PetscBool magnit,rel;
343: EPSWhich which;
345: PetscFunctionBegin;
346: *reason = EPS_CONVERGED_ITERATING;
347: firstev = ((EPSStoppingCtx)ctx)->firstev;
348: lastev = ((EPSStoppingCtx)ctx)->lastev;
349: thres = ((EPSStoppingCtx)ctx)->thres;
350: rel = ((EPSStoppingCtx)ctx)->threlative;
351: which = ((EPSStoppingCtx)ctx)->which;
352: magnit = (which==EPS_SMALLEST_MAGNITUDE || which==EPS_LARGEST_MAGNITUDE || which==EPS_TARGET_MAGNITUDE)? PETSC_TRUE: PETSC_FALSE;
353: if (nconv && magnit && which==EPS_TARGET_MAGNITUDE && ((rel && ((thres>1.0 && lastev>thres*firstev) || (thres<1.0 && lastev<thres*firstev))) || (!rel && lastev>thres))) {
354: if (!rel) PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: the eigenvalue magnitude %g is above the threshold %g\n",(double)lastev,(double)thres));
355: else if (thres>1.0) PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: the ratio %g/%g is above the threshold %g\n",(double)lastev,(double)firstev,(double)thres));
356: else PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: the ratio %g/%g is below the threshold %g\n",(double)lastev,(double)firstev,(double)thres));
357: *reason = EPS_CONVERGED_TOL;
358: } else if (nconv && magnit && ((which==EPS_LARGEST_MAGNITUDE && ((rel && lastev<thres*firstev) || (!rel && lastev<thres))) || (which==EPS_SMALLEST_MAGNITUDE && lastev>thres))) {
359: if (which==EPS_SMALLEST_MAGNITUDE) PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: the eigenvalue magnitude %g is above the threshold %g\n",(double)lastev,(double)thres));
360: else if (!rel) PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: the eigenvalue magnitude %g is below the threshold %g\n",(double)lastev,(double)thres));
361: else PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: the ratio %g/%g is below the threshold %g\n",(double)lastev,(double)firstev,(double)thres));
362: *reason = EPS_CONVERGED_TOL;
363: } else if (nconv && !magnit && ((which==EPS_LARGEST_REAL && lastev<thres) || (which==EPS_SMALLEST_REAL && lastev>thres))) {
364: if (which==EPS_LARGEST_REAL) PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: eigenvalue %g is below the threshold %g\n",(double)lastev,(double)thres));
365: else PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: eigenvalue %g is above the threshold %g\n",(double)lastev,(double)thres));
366: *reason = EPS_CONVERGED_TOL;
367: } else if (nev && nconv >= nev) {
368: PetscCall(PetscInfo(eps,"Linear eigensolver finished successfully: %" PetscInt_FMT " eigenpairs converged at iteration %" PetscInt_FMT "\n",nconv,its));
369: *reason = EPS_CONVERGED_TOL;
370: } else if (its >= max_it) {
371: *reason = EPS_DIVERGED_ITS;
372: PetscCall(PetscInfo(eps,"Linear eigensolver iteration reached maximum number of iterations (%" PetscInt_FMT ")\n",its));
373: }
374: PetscFunctionReturn(PETSC_SUCCESS);
375: }
377: /*
378: EPSComputeRitzVector - Computes the current Ritz vector.
380: Simple case (complex scalars or real scalars with Zi=NULL):
381: x = V*Zr (V is a basis of nv vectors, Zr has length nv)
383: Split case:
384: x = V*Zr y = V*Zi (Zr and Zi have length nv)
385: */
386: PetscErrorCode EPSComputeRitzVector(EPS eps,PetscScalar *Zr,PetscScalar *Zi,BV V,Vec x,Vec y)
387: {
388: PetscInt l,k;
389: PetscReal norm;
390: #if !defined(PETSC_USE_COMPLEX)
391: Vec z;
392: #endif
394: PetscFunctionBegin;
395: /* compute eigenvector */
396: PetscCall(BVGetActiveColumns(V,&l,&k));
397: PetscCall(BVSetActiveColumns(V,0,k));
398: PetscCall(BVMultVec(V,1.0,0.0,x,Zr));
400: /* purify eigenvector if necessary */
401: if (eps->purify) {
402: PetscCall(STApply(eps->st,x,y));
403: if (eps->ishermitian) PetscCall(BVNormVec(eps->V,y,NORM_2,&norm));
404: else PetscCall(VecNorm(y,NORM_2,&norm));
405: PetscCall(VecScale(y,1.0/norm));
406: PetscCall(VecCopy(y,x));
407: }
408: /* fix eigenvector if balancing is used */
409: if (!eps->ishermitian && eps->balance!=EPS_BALANCE_NONE && eps->D) PetscCall(VecPointwiseDivide(x,x,eps->D));
410: #if !defined(PETSC_USE_COMPLEX)
411: /* compute imaginary part of eigenvector */
412: if (Zi) {
413: PetscCall(BVMultVec(V,1.0,0.0,y,Zi));
414: if (eps->ispositive) {
415: PetscCall(BVCreateVec(V,&z));
416: PetscCall(STApply(eps->st,y,z));
417: PetscCall(VecNorm(z,NORM_2,&norm));
418: PetscCall(VecScale(z,1.0/norm));
419: PetscCall(VecCopy(z,y));
420: PetscCall(VecDestroy(&z));
421: }
422: if (eps->balance!=EPS_BALANCE_NONE && eps->D) PetscCall(VecPointwiseDivide(y,y,eps->D));
423: } else
424: #endif
425: PetscCall(VecSet(y,0.0));
427: /* normalize eigenvectors (when using balancing) */
428: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
429: #if !defined(PETSC_USE_COMPLEX)
430: if (Zi) PetscCall(VecNormalizeComplex(x,y,PETSC_TRUE,NULL));
431: else
432: #endif
433: PetscCall(VecNormalize(x,NULL));
434: }
435: PetscCall(BVSetActiveColumns(V,l,k));
436: PetscFunctionReturn(PETSC_SUCCESS);
437: }
439: /*
440: EPSBuildBalance_Krylov - uses a Krylov subspace method to compute the
441: diagonal matrix to be applied for balancing in non-Hermitian problems.
442: */
443: PetscErrorCode EPSBuildBalance_Krylov(EPS eps)
444: {
445: Vec z,p,r;
446: PetscInt i,j;
447: PetscReal norma;
448: PetscScalar *pz,*pD;
449: const PetscScalar *pr,*pp;
450: PetscRandom rand;
452: PetscFunctionBegin;
453: PetscCall(EPSSetWorkVecs(eps,3));
454: PetscCall(BVGetRandomContext(eps->V,&rand));
455: r = eps->work[0];
456: p = eps->work[1];
457: z = eps->work[2];
458: PetscCall(VecSet(eps->D,1.0));
460: for (j=0;j<eps->balance_its;j++) {
462: /* Build a random vector of +-1's */
463: PetscCall(VecSetRandom(z,rand));
464: PetscCall(VecGetArray(z,&pz));
465: for (i=0;i<eps->nloc;i++) {
466: if (PetscRealPart(pz[i])<0.5) pz[i]=-1.0;
467: else pz[i]=1.0;
468: }
469: PetscCall(VecRestoreArray(z,&pz));
471: /* Compute p=DA(D\z) */
472: PetscCall(VecPointwiseDivide(r,z,eps->D));
473: PetscCall(STApply(eps->st,r,p));
474: PetscCall(VecPointwiseMult(p,p,eps->D));
475: if (eps->balance == EPS_BALANCE_TWOSIDE) {
476: if (j==0) {
477: /* Estimate the matrix inf-norm */
478: PetscCall(VecAbs(p));
479: PetscCall(VecMax(p,NULL,&norma));
480: }
481: /* Compute r=D\(A'Dz) */
482: PetscCall(VecPointwiseMult(z,z,eps->D));
483: PetscCall(STApplyHermitianTranspose(eps->st,z,r));
484: PetscCall(VecPointwiseDivide(r,r,eps->D));
485: }
487: /* Adjust values of D */
488: PetscCall(VecGetArrayRead(r,&pr));
489: PetscCall(VecGetArrayRead(p,&pp));
490: PetscCall(VecGetArray(eps->D,&pD));
491: for (i=0;i<eps->nloc;i++) {
492: if (eps->balance == EPS_BALANCE_TWOSIDE) {
493: if (PetscAbsScalar(pp[i])>eps->balance_cutoff*norma && pr[i]!=0.0)
494: pD[i] *= PetscSqrtReal(PetscAbsScalar(pr[i]/pp[i]));
495: } else {
496: if (pp[i]!=0.0) pD[i] /= PetscAbsScalar(pp[i]);
497: }
498: }
499: PetscCall(VecRestoreArrayRead(r,&pr));
500: PetscCall(VecRestoreArrayRead(p,&pp));
501: PetscCall(VecRestoreArray(eps->D,&pD));
502: }
503: PetscFunctionReturn(PETSC_SUCCESS);
504: }