1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: EPS routines related to the solution process
12: */
14: #include <slepc/private/epsimpl.h> 15: #include <slepc/private/bvimpl.h> 16: #include <petscdraw.h>
18: PetscErrorCode EPSComputeVectors(EPS eps) 19: {
23: EPSCheckSolved(eps,1);
24: if (eps->state==EPS_STATE_SOLVED && eps->ops->computevectors) {
25: (*eps->ops->computevectors)(eps);
26: }
27: eps->state = EPS_STATE_EIGENVECTORS;
28: return(0);
29: }
31: #define SWAP(a,b,t) {t=a;a=b;b=t;} 33: static PetscErrorCode EPSComputeValues(EPS eps) 34: {
36: PetscBool injective,iscomp,isfilter;
37: PetscInt i,n,aux,nconv0;
38: Mat A,B=NULL,G,Z;
41: switch (eps->categ) {
42: case EPS_CATEGORY_KRYLOV:
43: case EPS_CATEGORY_OTHER:
44: STIsInjective(eps->st,&injective);
45: if (injective) {
46: /* one-to-one mapping: backtransform eigenvalues */
47: if (eps->ops->backtransform) {
48: (*eps->ops->backtransform)(eps);
49: } else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, spectral transform should have a backtransform operation");
50: } else {
51: /* compute eigenvalues from Rayleigh quotient */
52: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
53: if (!n) break;
54: EPSGetOperators(eps,&A,&B);
55: BVSetActiveColumns(eps->V,0,n);
56: DSGetCompact(eps->ds,&iscomp);
57: DSSetCompact(eps->ds,PETSC_FALSE);
58: DSGetMat(eps->ds,DS_MAT_A,&G);
59: BVMatProject(eps->V,A,eps->V,G);
60: DSRestoreMat(eps->ds,DS_MAT_A,&G);
61: if (B) {
62: DSGetMat(eps->ds,DS_MAT_B,&G);
63: BVMatProject(eps->V,B,eps->V,G);
64: DSRestoreMat(eps->ds,DS_MAT_A,&G);
65: }
66: DSSolve(eps->ds,eps->eigr,eps->eigi);
67: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
68: DSSynchronize(eps->ds,eps->eigr,eps->eigi);
69: DSSetCompact(eps->ds,iscomp);
70: if (eps->ishermitian && (!eps->isgeneralized || eps->ispositive)) { /* V = V * Z */
71: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
72: DSGetMat(eps->ds,DS_MAT_X,&Z);
73: BVMultInPlace(eps->V,Z,0,n);
74: MatDestroy(&Z);
75: }
76: /* in case of STFILTER discard computed eigenvalues that lie outside the wanted interval */
77: PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilter);
78: if (isfilter) {
79: nconv0 = eps->nconv;
80: for (i=0;i<eps->nconv;i++) {
81: if (PetscRealPart(eps->eigr[eps->perm[i]])<eps->inta || PetscRealPart(eps->eigr[eps->perm[i]])>eps->intb) {
82: eps->nconv--;
83: if (i<eps->nconv) { SWAP(eps->perm[i],eps->perm[eps->nconv],aux); i--; }
84: }
85: }
86: if (nconv0>eps->nconv) {
87: PetscInfo1(eps,"Discarded %D computed eigenvalues lying outside the interval\n",nconv0-eps->nconv);
88: }
89: }
90: }
91: break;
92: case EPS_CATEGORY_PRECOND:
93: case EPS_CATEGORY_CONTOUR:
94: /* eigenvalues already available as an output of the solver */
95: break;
96: }
97: return(0);
98: }
100: /*@
101: EPSSolve - Solves the eigensystem.
103: Collective on eps
105: Input Parameter:
106: . eps - eigensolver context obtained from EPSCreate()
108: Options Database Keys:
109: + -eps_view - print information about the solver used
110: . -eps_view_mat0 binary - save the first matrix (A) to the default binary viewer
111: . -eps_view_mat1 binary - save the second matrix (B) to the default binary viewer
112: . -eps_view_vectors binary - save the computed eigenvectors to the default binary viewer
113: . -eps_view_values - print computed eigenvalues
114: . -eps_converged_reason - print reason for convergence, and number of iterations
115: . -eps_error_absolute - print absolute errors of each eigenpair
116: . -eps_error_relative - print relative errors of each eigenpair
117: - -eps_error_backward - print backward errors of each eigenpair
119: Level: beginner
121: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
122: @*/
123: PetscErrorCode EPSSolve(EPS eps)124: {
126: PetscInt i;
127: STMatMode matmode;
128: Mat A,B;
132: if (eps->state>=EPS_STATE_SOLVED) return(0);
133: PetscLogEventBegin(EPS_Solve,eps,0,0,0);
135: /* call setup */
136: EPSSetUp(eps);
137: eps->nconv = 0;
138: eps->its = 0;
139: for (i=0;i<eps->ncv;i++) {
140: eps->eigr[i] = 0.0;
141: eps->eigi[i] = 0.0;
142: eps->errest[i] = 0.0;
143: eps->perm[i] = i;
144: }
145: EPSViewFromOptions(eps,NULL,"-eps_view_pre");
146: RGViewFromOptions(eps->rg,NULL,"-rg_view");
148: /* call solver */
149: (*eps->ops->solve)(eps);
150: eps->state = EPS_STATE_SOLVED;
152: STGetMatMode(eps->st,&matmode);
153: if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
154: /* Purify eigenvectors before reverting operator */
155: EPSComputeVectors(eps);
156: }
157: STPostSolve(eps->st);
159: if (!eps->reason) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
161: /* Map eigenvalues back to the original problem if appropriate */
162: EPSComputeValues(eps);
164: #if !defined(PETSC_USE_COMPLEX)
165: /* reorder conjugate eigenvalues (positive imaginary first) */
166: for (i=0;i<eps->nconv-1;i++) {
167: if (eps->eigi[i] != 0) {
168: if (eps->eigi[i] < 0) {
169: eps->eigi[i] = -eps->eigi[i];
170: eps->eigi[i+1] = -eps->eigi[i+1];
171: /* the next correction only works with eigenvectors */
172: EPSComputeVectors(eps);
173: BVScaleColumn(eps->V,i+1,-1.0);
174: }
175: i++;
176: }
177: }
178: #endif
180: /* sort eigenvalues according to eps->which parameter */
181: SlepcSortEigenvalues(eps->sc,eps->nconv,eps->eigr,eps->eigi,eps->perm);
182: PetscLogEventEnd(EPS_Solve,eps,0,0,0);
184: /* various viewers */
185: EPSViewFromOptions(eps,NULL,"-eps_view");
186: EPSReasonViewFromOptions(eps);
187: EPSErrorViewFromOptions(eps);
188: EPSValuesViewFromOptions(eps);
189: EPSVectorsViewFromOptions(eps);
190: EPSGetOperators(eps,&A,&B);
191: MatViewFromOptions(A,(PetscObject)eps,"-eps_view_mat0");
192: if (eps->isgeneralized) {
193: MatViewFromOptions(B,(PetscObject)eps,"-eps_view_mat1");
194: }
196: /* Remove deflation and initial subspaces */
197: if (eps->nds) {
198: BVSetNumConstraints(eps->V,0);
199: eps->nds = 0;
200: }
201: eps->nini = 0;
202: return(0);
203: }
205: /*@
206: EPSGetIterationNumber - Gets the current iteration number. If the
207: call to EPSSolve() is complete, then it returns the number of iterations
208: carried out by the solution method.
210: Not Collective
212: Input Parameter:
213: . eps - the eigensolver context
215: Output Parameter:
216: . its - number of iterations
218: Note:
219: During the i-th iteration this call returns i-1. If EPSSolve() is
220: complete, then parameter "its" contains either the iteration number at
221: which convergence was successfully reached, or failure was detected.
222: Call EPSGetConvergedReason() to determine if the solver converged or
223: failed and why.
225: Level: intermediate
227: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
228: @*/
229: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)230: {
234: *its = eps->its;
235: return(0);
236: }
238: /*@
239: EPSGetConverged - Gets the number of converged eigenpairs.
241: Not Collective
243: Input Parameter:
244: . eps - the eigensolver context
246: Output Parameter:
247: . nconv - number of converged eigenpairs
249: Note:
250: This function should be called after EPSSolve() has finished.
252: Level: beginner
254: .seealso: EPSSetDimensions(), EPSSolve()
255: @*/
256: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)257: {
261: EPSCheckSolved(eps,1);
262: *nconv = eps->nconv;
263: return(0);
264: }
266: /*@
267: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
268: stopped.
270: Not Collective
272: Input Parameter:
273: . eps - the eigensolver context
275: Output Parameter:
276: . reason - negative value indicates diverged, positive value converged
278: Notes:
279: Possible values for reason are
280: + EPS_CONVERGED_TOL - converged up to tolerance
281: . EPS_CONVERGED_USER - converged due to a user-defined condition
282: . EPS_DIVERGED_ITS - required more than max_it iterations to reach convergence
283: . EPS_DIVERGED_BREAKDOWN - generic breakdown in method
284: - EPS_DIVERGED_SYMMETRY_LOST - pseudo-Lanczos was not able to keep symmetry
286: Can only be called after the call to EPSSolve() is complete.
288: Level: intermediate
290: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason291: @*/
292: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)293: {
297: EPSCheckSolved(eps,1);
298: *reason = eps->reason;
299: return(0);
300: }
302: /*@
303: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
304: subspace.
306: Not Collective, but vectors are shared by all processors that share the EPS308: Input Parameter:
309: . eps - the eigensolver context
311: Output Parameter:
312: . v - an array of vectors
314: Notes:
315: This function should be called after EPSSolve() has finished.
317: The user should provide in v an array of nconv vectors, where nconv is
318: the value returned by EPSGetConverged().
320: The first k vectors returned in v span an invariant subspace associated
321: with the first k computed eigenvalues (note that this is not true if the
322: k-th eigenvalue is complex and matrix A is real; in this case the first
323: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
324: in X for all x in X (a similar definition applies for generalized
325: eigenproblems).
327: Level: intermediate
329: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve()
330: @*/
331: PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec v[])332: {
334: PetscInt i;
335: BV V=eps->V;
336: Vec w;
342: EPSCheckSolved(eps,1);
343: if (!eps->ishermitian && eps->state==EPS_STATE_EIGENVECTORS) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector or EPSComputeError");
344: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
345: BVDuplicateResize(eps->V,eps->nconv,&V);
346: BVSetActiveColumns(eps->V,0,eps->nconv);
347: BVCopy(eps->V,V);
348: for (i=0;i<eps->nconv;i++) {
349: BVGetColumn(V,i,&w);
350: VecPointwiseDivide(w,w,eps->D);
351: BVRestoreColumn(V,i,&w);
352: }
353: BVOrthogonalize(V,NULL);
354: }
355: for (i=0;i<eps->nconv;i++) {
356: BVCopyVec(V,i,v[i]);
357: }
358: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
359: BVDestroy(&V);
360: }
361: return(0);
362: }
364: /*@C
365: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
366: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
368: Logically Collective on eps
370: Input Parameters:
371: + eps - eigensolver context
372: - i - index of the solution
374: Output Parameters:
375: + eigr - real part of eigenvalue
376: . eigi - imaginary part of eigenvalue
377: . Vr - real part of eigenvector
378: - Vi - imaginary part of eigenvector
380: Notes:
381: It is allowed to pass NULL for Vr and Vi, if the eigenvector is not
382: required. Otherwise, the caller must provide valid Vec objects, i.e.,
383: they must be created by the calling program with e.g. MatCreateVecs().
385: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
386: configured with complex scalars the eigenvalue is stored
387: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
388: set to zero). In both cases, the user can pass NULL in eigi and Vi.
390: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
391: Eigenpairs are indexed according to the ordering criterion established
392: with EPSSetWhichEigenpairs().
394: The 2-norm of the eigenvector is one unless the problem is generalized
395: Hermitian. In this case the eigenvector is normalized with respect to the
396: norm defined by the B matrix.
398: Level: beginner
400: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSGetLeftEigenvector(), EPSSolve(),
401: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
402: @*/
403: PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)404: {
410: EPSCheckSolved(eps,1);
411: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
412: EPSGetEigenvalue(eps,i,eigr,eigi);
413: if (Vr || Vi) { EPSGetEigenvector(eps,i,Vr,Vi); }
414: return(0);
415: }
417: /*@C
418: EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve().
420: Not Collective
422: Input Parameters:
423: + eps - eigensolver context
424: - i - index of the solution
426: Output Parameters:
427: + eigr - real part of eigenvalue
428: - eigi - imaginary part of eigenvalue
430: Notes:
431: If the eigenvalue is real, then eigi is set to zero. If PETSc is
432: configured with complex scalars the eigenvalue is stored
433: directly in eigr (eigi is set to zero).
435: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
436: Eigenpairs are indexed according to the ordering criterion established
437: with EPSSetWhichEigenpairs().
439: Level: beginner
441: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
442: @*/
443: PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi)444: {
445: PetscInt k;
449: EPSCheckSolved(eps,1);
450: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
451: k = eps->perm[i];
452: #if defined(PETSC_USE_COMPLEX)
453: if (eigr) *eigr = eps->eigr[k];
454: if (eigi) *eigi = 0;
455: #else
456: if (eigr) *eigr = eps->eigr[k];
457: if (eigi) *eigi = eps->eigi[k];
458: #endif
459: return(0);
460: }
462: /*@
463: EPSGetEigenvector - Gets the i-th right eigenvector as computed by EPSSolve().
465: Logically Collective on eps
467: Input Parameters:
468: + eps - eigensolver context
469: - i - index of the solution
471: Output Parameters:
472: + Vr - real part of eigenvector
473: - Vi - imaginary part of eigenvector
475: Notes:
476: The caller must provide valid Vec objects, i.e., they must be created
477: by the calling program with e.g. MatCreateVecs().
479: If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
480: configured with complex scalars the eigenvector is stored
481: directly in Vr (Vi is set to zero). In any case, the user can pass NULL in Vr
482: or Vi if one of them is not required.
484: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
485: Eigenpairs are indexed according to the ordering criterion established
486: with EPSSetWhichEigenpairs().
488: The 2-norm of the eigenvector is one unless the problem is generalized
489: Hermitian. In this case the eigenvector is normalized with respect to the
490: norm defined by the B matrix.
492: Level: beginner
494: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair(), EPSGetLeftEigenvector()
495: @*/
496: PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi)497: {
499: PetscInt k;
506: EPSCheckSolved(eps,1);
507: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
508: EPSComputeVectors(eps);
509: k = eps->perm[i];
510: BV_GetEigenvector(eps->V,k,eps->eigi[k],Vr,Vi);
511: return(0);
512: }
514: /*@
515: EPSGetLeftEigenvector - Gets the i-th left eigenvector as computed by EPSSolve().
517: Logically Collective on eps
519: Input Parameters:
520: + eps - eigensolver context
521: - i - index of the solution
523: Output Parameters:
524: + Wr - real part of left eigenvector
525: - Wi - imaginary part of left eigenvector
527: Notes:
528: The caller must provide valid Vec objects, i.e., they must be created
529: by the calling program with e.g. MatCreateVecs().
531: If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is
532: configured with complex scalars the eigenvector is stored directly in Wr
533: (Wi is set to zero). In any case, the user can pass NULL in Wr or Wi if
534: one of them is not required.
536: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
537: Eigensolutions are indexed according to the ordering criterion established
538: with EPSSetWhichEigenpairs().
540: Left eigenvectors are available only if the twosided flag was set, see
541: EPSSetTwoSided().
543: Level: intermediate
545: .seealso: EPSGetEigenvector(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSSetTwoSided()
546: @*/
547: PetscErrorCode EPSGetLeftEigenvector(EPS eps,PetscInt i,Vec Wr,Vec Wi)548: {
550: PetscInt k;
557: EPSCheckSolved(eps,1);
558: if (!eps->twosided) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with EPSSetTwoSided");
559: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
560: EPSComputeVectors(eps);
561: k = eps->perm[i];
562: BV_GetEigenvector(eps->W,k,eps->eigi[k],Wr,Wi);
563: return(0);
564: }
566: /*@
567: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
568: computed eigenpair.
570: Not Collective
572: Input Parameter:
573: + eps - eigensolver context
574: - i - index of eigenpair
576: Output Parameter:
577: . errest - the error estimate
579: Notes:
580: This is the error estimate used internally by the eigensolver. The actual
581: error bound can be computed with EPSComputeError(). See also the users
582: manual for details.
584: Level: advanced
586: .seealso: EPSComputeError()
587: @*/
588: PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)589: {
593: EPSCheckSolved(eps,1);
594: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
595: *errest = eps->errest[eps->perm[i]];
596: return(0);
597: }
599: /*
600: EPSComputeResidualNorm_Private - Computes the norm of the residual vector
601: associated with an eigenpair.
603: Input Parameters:
604: trans - whether A' must be used instead of A
605: kr,ki - eigenvalue
606: xr,xi - eigenvector
607: z - three work vectors (the second one not referenced in complex scalars)
608: */
609: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)610: {
612: PetscInt nmat;
613: Mat A,B;
614: Vec u,w;
615: PetscScalar alpha;
616: #if !defined(PETSC_USE_COMPLEX)
617: Vec v;
618: PetscReal ni,nr;
619: #endif
620: PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultHermitianTranspose: MatMult;
623: u = z[0]; w = z[2];
624: STGetNumMatrices(eps->st,&nmat);
625: STGetMatrix(eps->st,0,&A);
626: if (nmat>1) { STGetMatrix(eps->st,1,&B); }
628: #if !defined(PETSC_USE_COMPLEX)
629: v = z[1];
630: if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
631: #endif
632: (*matmult)(A,xr,u); /* u=A*x */
633: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
634: if (nmat>1) { (*matmult)(B,xr,w); }
635: else { VecCopy(xr,w); } /* w=B*x */
636: alpha = trans? -PetscConj(kr): -kr;
637: VecAXPY(u,alpha,w); /* u=A*x-k*B*x */
638: }
639: VecNorm(u,NORM_2,norm);
640: #if !defined(PETSC_USE_COMPLEX)
641: } else {
642: (*matmult)(A,xr,u); /* u=A*xr */
643: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
644: if (nmat>1) { (*matmult)(B,xr,v); }
645: else { VecCopy(xr,v); } /* v=B*xr */
646: VecAXPY(u,-kr,v); /* u=A*xr-kr*B*xr */
647: if (nmat>1) { (*matmult)(B,xi,w); }
648: else { VecCopy(xi,w); } /* w=B*xi */
649: VecAXPY(u,trans?-ki:ki,w); /* u=A*xr-kr*B*xr+ki*B*xi */
650: }
651: VecNorm(u,NORM_2,&nr);
652: (*matmult)(A,xi,u); /* u=A*xi */
653: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
654: VecAXPY(u,-kr,w); /* u=A*xi-kr*B*xi */
655: VecAXPY(u,trans?ki:-ki,v); /* u=A*xi-kr*B*xi-ki*B*xr */
656: }
657: VecNorm(u,NORM_2,&ni);
658: *norm = SlepcAbsEigenvalue(nr,ni);
659: }
660: #endif
661: return(0);
662: }
664: /*@
665: EPSComputeError - Computes the error (based on the residual norm) associated
666: with the i-th computed eigenpair.
668: Collective on eps
670: Input Parameter:
671: + eps - the eigensolver context
672: . i - the solution index
673: - type - the type of error to compute
675: Output Parameter:
676: . error - the error
678: Notes:
679: The error can be computed in various ways, all of them based on the residual
680: norm ||Ax-kBx||_2 where k is the eigenvalue and x is the eigenvector.
682: Level: beginner
684: .seealso: EPSErrorType, EPSSolve(), EPSGetErrorEstimate()
685: @*/
686: PetscErrorCode EPSComputeError(EPS eps,PetscInt i,EPSErrorType type,PetscReal *error)687: {
689: Mat A,B;
690: Vec xr,xi,w[3];
691: PetscReal t,vecnorm=1.0,errorl;
692: PetscScalar kr,ki;
693: PetscBool flg;
700: EPSCheckSolved(eps,1);
702: /* allocate work vectors */
703: #if defined(PETSC_USE_COMPLEX)
704: EPSSetWorkVecs(eps,3);
705: xi = NULL;
706: w[1] = NULL;
707: #else
708: EPSSetWorkVecs(eps,5);
709: xi = eps->work[3];
710: w[1] = eps->work[4];
711: #endif
712: xr = eps->work[0];
713: w[0] = eps->work[1];
714: w[2] = eps->work[2];
716: /* compute residual norm */
717: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
718: EPSComputeResidualNorm_Private(eps,PETSC_FALSE,kr,ki,xr,xi,w,error);
720: /* compute 2-norm of eigenvector */
721: if (eps->problem_type==EPS_GHEP) {
722: VecNorm(xr,NORM_2,&vecnorm);
723: }
725: /* if two-sided, compute left residual norm and take the maximum */
726: if (eps->twosided) {
727: EPSGetLeftEigenvector(eps,i,xr,xi);
728: EPSComputeResidualNorm_Private(eps,PETSC_TRUE,kr,ki,xr,xi,w,&errorl);
729: *error = PetscMax(*error,errorl);
730: }
732: /* compute error */
733: switch (type) {
734: case EPS_ERROR_ABSOLUTE:
735: break;
736: case EPS_ERROR_RELATIVE:
737: *error /= SlepcAbsEigenvalue(kr,ki)*vecnorm;
738: break;
739: case EPS_ERROR_BACKWARD:
740: /* initialization of matrix norms */
741: if (!eps->nrma) {
742: STGetMatrix(eps->st,0,&A);
743: MatHasOperation(A,MATOP_NORM,&flg);
744: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
745: MatNorm(A,NORM_INFINITY,&eps->nrma);
746: }
747: if (eps->isgeneralized) {
748: if (!eps->nrmb) {
749: STGetMatrix(eps->st,1,&B);
750: MatHasOperation(B,MATOP_NORM,&flg);
751: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
752: MatNorm(B,NORM_INFINITY,&eps->nrmb);
753: }
754: } else eps->nrmb = 1.0;
755: t = SlepcAbsEigenvalue(kr,ki);
756: *error /= (eps->nrma+t*eps->nrmb)*vecnorm;
757: break;
758: default:759: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
760: }
761: return(0);
762: }
764: /*
765: EPSGetStartVector - Generate a suitable vector to be used as the starting vector
766: for the recurrence that builds the right subspace.
768: Collective on eps
770: Input Parameters:
771: + eps - the eigensolver context
772: - i - iteration number
774: Output Parameters:
775: . breakdown - flag indicating that a breakdown has occurred
777: Notes:
778: The start vector is computed from another vector: for the first step (i=0),
779: the first initial vector is used (see EPSSetInitialSpace()); otherwise a random
780: vector is created. Then this vector is forced to be in the range of OP (only
781: for generalized definite problems) and orthonormalized with respect to all
782: V-vectors up to i-1. The resulting vector is placed in V[i].
784: The flag breakdown is set to true if either i=0 and the vector belongs to the
785: deflation space, or i>0 and the vector is linearly dependent with respect
786: to the V-vectors.
787: */
788: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,PetscBool *breakdown)789: {
791: PetscReal norm;
792: PetscBool lindep;
793: Vec w,z;
799: /* For the first step, use the first initial vector, otherwise a random one */
800: if (i>0 || eps->nini==0) {
801: BVSetRandomColumn(eps->V,i);
802: }
804: /* Force the vector to be in the range of OP for definite generalized problems */
805: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
806: BVCreateVec(eps->V,&w);
807: BVCopyVec(eps->V,i,w);
808: BVGetColumn(eps->V,i,&z);
809: STApply(eps->st,w,z);
810: BVRestoreColumn(eps->V,i,&z);
811: VecDestroy(&w);
812: }
814: /* Orthonormalize the vector with respect to previous vectors */
815: BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&lindep);
816: if (breakdown) *breakdown = lindep;
817: else if (lindep || norm == 0.0) {
818: if (i==0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Initial vector is zero or belongs to the deflation space");
819: else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to generate more start vectors");
820: }
821: BVScaleColumn(eps->V,i,1.0/norm);
822: return(0);
823: }
825: /*
826: EPSGetLeftStartVector - Generate a suitable vector to be used as the left starting
827: vector for the recurrence that builds the left subspace. See EPSGetStartVector().
828: */
829: PetscErrorCode EPSGetLeftStartVector(EPS eps,PetscInt i,PetscBool *breakdown)830: {
832: PetscReal norm;
833: PetscBool lindep;
839: /* For the first step, use the first initial vector, otherwise a random one */
840: if (i>0 || eps->ninil==0) {
841: BVSetRandomColumn(eps->W,i);
842: }
844: /* Orthonormalize the vector with respect to previous vectors */
845: BVOrthogonalizeColumn(eps->W,i,NULL,&norm,&lindep);
846: if (breakdown) *breakdown = lindep;
847: else if (lindep || norm == 0.0) {
848: if (i==0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Left initial vector is zero");
849: else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to generate more left start vectors");
850: }
851: BVScaleColumn(eps->W,i,1.0/norm);
852: return(0);
853: }