Actual source code: solve.c
1: /*
2: EPS routines related to the solution process.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2009, Universidad 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 private/epsimpl.h
28: /*@
29: EPSSolve - Solves the eigensystem.
31: Collective on EPS
33: Input Parameter:
34: . eps - eigensolver context obtained from EPSCreate()
36: Options Database:
37: + -eps_view - print information about the solver used
38: . -eps_view_binary - save the matrices to the default binary file
39: - -eps_plot_eigs - plot computed eigenvalues
41: Level: beginner
43: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
44: @*/
45: PetscErrorCode EPSSolve(EPS eps)
46: {
48: PetscInt i;
49: PetscReal re,im;
50: PetscTruth flg;
51: PetscViewer viewer;
52: PetscDraw draw;
53: PetscDrawSP drawsp;
54: STMatMode matmode;
59: PetscOptionsHasName(((PetscObject)eps)->prefix,"-eps_view_binary",&flg);
60: if (flg) {
61: Mat A,B;
62: STGetOperators(eps->OP,&A,&B);
63: MatView(A,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
64: if (B) MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
65: }
67: /* reset the convergence flag from the previous solves */
68: eps->reason = EPS_CONVERGED_ITERATING;
70: if (!eps->setupcalled){ EPSSetUp(eps); }
71: STResetOperationCounters(eps->OP);
72: IPResetOperationCounters(eps->ip);
73: eps->evecsavailable = PETSC_FALSE;
74: eps->nconv = 0;
75: eps->its = 0;
76: for (i=0;i<eps->ncv;i++) eps->eigr[i]=eps->eigi[i]=eps->errest[i]=0.0;
77: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->ncv);
79: PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
81: switch (eps->solverclass) {
82: case EPS_ONE_SIDE:
83: (*eps->ops->solve)(eps); break;
84: case EPS_TWO_SIDE:
85: if (eps->ops->solvets) {
86: (*eps->ops->solvets)(eps); break;
87: } else {
88: SETERRQ(1,"Two-sided version unavailable for this solver");
89: }
90: default:
91: SETERRQ(1,"Wrong value of eps->solverclass");
92: }
94: STGetMatMode(eps->OP,&matmode);
95: if (matmode == STMATMODE_INPLACE && eps->ispositive) {
96: /* Purge eigenvectors before reverting operator */
97: (*eps->ops->computevectors)(eps);
98: }
99: STPostSolve(eps->OP);
100: PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
102: if (!eps->reason) {
103: SETERRQ(1,"Internal error, solver returned without setting converged reason");
104: }
106: /* Map eigenvalues back to the original problem, necessary in some
107: * spectral transformations */
108: (*eps->ops->backtransform)(eps);
110: /* Adjust left eigenvectors in generalized problems: y = B^T y */
111: if (eps->isgeneralized && eps->solverclass == EPS_TWO_SIDE) {
112: Mat B;
113: KSP ksp;
114: Vec w;
115: STGetOperators(eps->OP,PETSC_NULL,&B);
116: KSPCreate(((PetscObject)eps)->comm,&ksp);
117: KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
118: KSPSetFromOptions(ksp);
119: MatGetVecs(B,PETSC_NULL,&w);
120: for (i=0;i<eps->nconv;i++) {
121: VecCopy(eps->W[i],w);
122: KSPSolveTranspose(ksp,w,eps->W[i]);
123: }
124: KSPDestroy(ksp);
125: VecDestroy(w);
126: }
127:
129: #ifndef PETSC_USE_COMPLEX
130: /* reorder conjugate eigenvalues (positive imaginary first) */
131: for (i=0; i<eps->nconv-1; i++) {
132: if (eps->eigi[i] != 0) {
133: if (eps->eigi[i] < 0) {
134: eps->eigi[i] = -eps->eigi[i];
135: eps->eigi[i+1] = -eps->eigi[i+1];
136: VecScale(eps->V[i+1],-1.0);
137: }
138: i++;
139: }
140: }
141: #endif
143: /* sort eigenvalues according to eps->which parameter */
144: PetscFree(eps->perm);
145: if (eps->nconv > 0) {
146: PetscMalloc(sizeof(PetscInt)*eps->nconv, &eps->perm);
147: EPSSortEigenvalues(eps->nconv, eps->eigr, eps->eigi, eps->which, eps->nconv, eps->perm);
148: }
150: PetscOptionsHasName(((PetscObject)eps)->prefix,"-eps_view",&flg);
151: if (flg && !PetscPreLoadingOn) { EPSView(eps,PETSC_VIEWER_STDOUT_WORLD); }
153: PetscOptionsHasName(((PetscObject)eps)->prefix,"-eps_plot_eigs",&flg);
154: if (flg) {
155: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
156: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
157: PetscViewerDrawGetDraw(viewer,0,&draw);
158: PetscDrawSPCreate(draw,1,&drawsp);
159: for( i=0; i<eps->nconv; i++ ) {
160: #if defined(PETSC_USE_COMPLEX)
161: re = PetscRealPart(eps->eigr[i]);
162: im = PetscImaginaryPart(eps->eigi[i]);
163: #else
164: re = eps->eigr[i];
165: im = eps->eigi[i];
166: #endif
167: PetscDrawSPAddPoint(drawsp,&re,&im);
168: }
169: PetscDrawSPDraw(drawsp);
170: PetscDrawSPDestroy(drawsp);
171: PetscViewerDestroy(viewer);
172: }
174: return(0);
175: }
179: /*@
180: EPSGetIterationNumber - Gets the current iteration number. If the
181: call to EPSSolve() is complete, then it returns the number of iterations
182: carried out by the solution method.
183:
184: Not Collective
186: Input Parameter:
187: . eps - the eigensolver context
189: Output Parameter:
190: . its - number of iterations
192: Level: intermediate
194: Note:
195: During the i-th iteration this call returns i-1. If EPSSolve() is
196: complete, then parameter "its" contains either the iteration number at
197: which convergence was successfully reached, or failure was detected.
198: Call EPSGetConvergedReason() to determine if the solver converged or
199: failed and why.
201: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
202: @*/
203: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
204: {
208: *its = eps->its;
209: return(0);
210: }
214: /*@
215: EPSGetOperationCounters - Gets the total number of operator applications,
216: inner product operations and linear iterations used by the ST object
217: during the last EPSSolve() call.
219: Not Collective
221: Input Parameter:
222: . eps - EPS context
224: Output Parameter:
225: + ops - number of operator applications
226: . dots - number of inner product operations
227: - lits - number of linear iterations
229: Notes:
230: When the eigensolver algorithm invokes STApply() then a linear system
231: must be solved (except in the case of standard eigenproblems and shift
232: transformation). The number of iterations required in this solve is
233: accumulated into a counter whose value is returned by this function.
235: These counters are reset to zero at each successive call to EPSSolve().
237: Level: intermediate
239: @*/
240: PetscErrorCode EPSGetOperationCounters(EPS eps,PetscInt* ops,PetscInt* dots,PetscInt* lits)
241: {
246: STGetOperationCounters(eps->OP,ops,lits);
247: if (dots) {
248: IPGetOperationCounters(eps->ip,dots);
249: }
250: return(0);
251: }
255: /*@
256: EPSGetConverged - Gets the number of converged eigenpairs.
258: Not Collective
260: Input Parameter:
261: . eps - the eigensolver context
262:
263: Output Parameter:
264: . nconv - number of converged eigenpairs
266: Note:
267: This function should be called after EPSSolve() has finished.
269: Level: beginner
271: .seealso: EPSSetDimensions()
272: @*/
273: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
274: {
278: *nconv = eps->nconv;
279: return(0);
280: }
285: /*@C
286: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
287: stopped.
289: Not Collective
291: Input Parameter:
292: . eps - the eigensolver context
294: Output Parameter:
295: . reason - negative value indicates diverged, positive value converged
296: (see EPSConvergedReason)
298: Possible values for reason:
299: + EPS_CONVERGED_TOL - converged up to tolerance
300: . EPS_DIVERGED_ITS - required more than its to reach convergence
301: . EPS_DIVERGED_BREAKDOWN - generic breakdown in method
302: - EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric
304: Level: intermediate
306: Notes: Can only be called after the call to EPSSolve() is complete.
308: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
309: @*/
310: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
311: {
315: *reason = eps->reason;
316: return(0);
317: }
321: /*@
322: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
323: subspace.
325: Not Collective
327: Input Parameter:
328: . eps - the eigensolver context
329:
330: Output Parameter:
331: . v - an array of vectors
333: Notes:
334: This function should be called after EPSSolve() has finished.
336: The user should provide in v an array of nconv vectors, where nconv is
337: the value returned by EPSGetConverged().
339: The first k vectors returned in v span an invariant subspace associated
340: with the first k computed eigenvalues (note that this is not true if the
341: k-th eigenvalue is complex and matrix A is real; in this case the first
342: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
343: in X for all x in X (a similar definition applies for generalized
344: eigenproblems).
346: Level: intermediate
348: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetLeftInvariantSubspace()
349: @*/
350: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
351: {
353: PetscInt i;
359: if (!eps->V) {
360: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
361: }
362: if (!eps->ishermitian && eps->evecsavailable) {
363: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetRightVector,EPSComputeRelativeError or EPSComputeResidualNorm");
364: }
365: for (i=0;i<eps->nconv;i++) {
366: VecCopy(eps->V[i],v[i]);
367: }
368: return(0);
369: }
373: /*@
374: EPSGetLeftInvariantSubspace - Gets an orthonormal basis of the computed left
375: invariant subspace (only available in two-sided eigensolvers).
377: Not Collective
379: Input Parameter:
380: . eps - the eigensolver context
381:
382: Output Parameter:
383: . v - an array of vectors
385: Notes:
386: This function should be called after EPSSolve() has finished.
388: The user should provide in v an array of nconv vectors, where nconv is
389: the value returned by EPSGetConverged().
391: The first k vectors returned in v span a left invariant subspace associated
392: with the first k computed eigenvalues (note that this is not true if the
393: k-th eigenvalue is complex and matrix A is real; in this case the first
394: k+1 vectors should be used). A left invariant subspace Y of A satisfies y'A
395: in Y for all y in Y (a similar definition applies for generalized
396: eigenproblems).
398: Level: intermediate
400: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
401: @*/
402: PetscErrorCode EPSGetLeftInvariantSubspace(EPS eps, Vec *v)
403: {
405: PetscInt i;
411: if (!eps->W) {
412: if (eps->solverclass!=EPS_TWO_SIDE) {
413: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
414: } else {
415: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
416: }
417: }
418: for (i=0;i<eps->nconv;i++) {
419: VecCopy(eps->W[i],v[i]);
420: }
421: return(0);
422: }
426: /*@
427: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
428: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
430: Not Collective
432: Input Parameters:
433: + eps - eigensolver context
434: - i - index of the solution
436: Output Parameters:
437: + eigr - real part of eigenvalue
438: . eigi - imaginary part of eigenvalue
439: . Vr - real part of eigenvector
440: - Vi - imaginary part of eigenvector
442: Notes:
443: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
444: configured with complex scalars the eigenvalue is stored
445: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
446: set to zero).
448: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
449: Eigenpairs are indexed according to the ordering criterion established
450: with EPSSetWhichEigenpairs().
452: Level: beginner
454: .seealso: EPSGetValue(), EPSGetRightVector(), EPSGetLeftVector(), EPSSolve(),
455: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
456: @*/
457: PetscErrorCode EPSGetEigenpair(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
458: {
463: if (!eps->eigr || !eps->eigi || !eps->V) {
464: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
465: }
466: if (i<0 || i>=eps->nconv) {
467: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
468: }
469: EPSGetValue(eps,i,eigr,eigi);
470: EPSGetRightVector(eps,i,Vr,Vi);
471:
472: return(0);
473: }
477: /*@
478: EPSGetValue - Gets the i-th eigenvalue as computed by EPSSolve().
480: Not Collective
482: Input Parameters:
483: + eps - eigensolver context
484: - i - index of the solution
486: Output Parameters:
487: + eigr - real part of eigenvalue
488: - eigi - imaginary part of eigenvalue
490: Notes:
491: If the eigenvalue is real, then eigi is set to zero. If PETSc is
492: configured with complex scalars the eigenvalue is stored
493: directly in eigr (eigi is set to zero).
495: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
496: Eigenpairs are indexed according to the ordering criterion established
497: with EPSSetWhichEigenpairs().
499: Level: beginner
501: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
502: EPSGetEigenpair()
503: @*/
504: PetscErrorCode EPSGetValue(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi)
505: {
506: PetscInt k;
510: if (!eps->eigr || !eps->eigi) {
511: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
512: }
513: if (i<0 || i>=eps->nconv) {
514: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
515: }
517: if (!eps->perm) k = i;
518: else k = eps->perm[i];
519: #ifdef PETSC_USE_COMPLEX
520: if (eigr) *eigr = eps->eigr[k];
521: if (eigi) *eigi = 0;
522: #else
523: if (eigr) *eigr = eps->eigr[k];
524: if (eigi) *eigi = eps->eigi[k];
525: #endif
526:
527: return(0);
528: }
532: /*@
533: EPSGetRightVector - Gets the i-th right eigenvector as computed by EPSSolve().
535: Not Collective
537: Input Parameters:
538: + eps - eigensolver context
539: - i - index of the solution
541: Output Parameters:
542: + Vr - real part of eigenvector
543: - Vi - imaginary part of eigenvector
545: Notes:
546: If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
547: configured with complex scalars the eigenvector is stored
548: directly in Vr (Vi is set to zero).
550: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
551: Eigenpairs are indexed according to the ordering criterion established
552: with EPSSetWhichEigenpairs().
554: Level: beginner
556: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
557: EPSGetEigenpair(), EPSGetLeftVector()
558: @*/
559: PetscErrorCode EPSGetRightVector(EPS eps, PetscInt i, Vec Vr, Vec Vi)
560: {
562: PetscInt k;
566: if (!eps->V) {
567: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
568: }
569: if (i<0 || i>=eps->nconv) {
570: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
571: }
572: if (!eps->evecsavailable && (Vr || Vi) ) {
573: (*eps->ops->computevectors)(eps);
574: }
576: if (!eps->perm) k = i;
577: else k = eps->perm[i];
578: #ifdef PETSC_USE_COMPLEX
579: if (Vr) { VecCopy(eps->V[k], Vr); }
580: if (Vi) { VecSet(Vi,0.0); }
581: #else
582: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
583: if (Vr) { VecCopy(eps->V[k], Vr); }
584: if (Vi) { VecCopy(eps->V[k+1], Vi); }
585: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
586: if (Vr) { VecCopy(eps->V[k-1], Vr); }
587: if (Vi) {
588: VecCopy(eps->V[k], Vi);
589: VecScale(Vi,-1.0);
590: }
591: } else { /* real eigenvalue */
592: if (Vr) { VecCopy(eps->V[k], Vr); }
593: if (Vi) { VecSet(Vi,0.0); }
594: }
595: #endif
596:
597: return(0);
598: }
602: /*@
603: EPSGetLeftVector - Gets the i-th left eigenvector as computed by EPSSolve()
604: (only available in two-sided eigensolvers).
606: Not Collective
608: Input Parameters:
609: + eps - eigensolver context
610: - i - index of the solution
612: Output Parameters:
613: + Wr - real part of eigenvector
614: - Wi - imaginary part of eigenvector
616: Notes:
617: If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is
618: configured with complex scalars the eigenvector is stored
619: directly in Wr (Wi is set to zero).
621: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
622: Eigenpairs are indexed according to the ordering criterion established
623: with EPSSetWhichEigenpairs().
625: Level: beginner
627: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
628: EPSGetEigenpair(), EPSGetLeftVector()
629: @*/
630: PetscErrorCode EPSGetLeftVector(EPS eps, PetscInt i, Vec Wr, Vec Wi)
631: {
633: PetscInt k;
637: if (!eps->W) {
638: if (eps->solverclass!=EPS_TWO_SIDE) {
639: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
640: } else {
641: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
642: }
643: }
644: if (i<0 || i>=eps->nconv) {
645: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
646: }
647: if (!eps->evecsavailable && (Wr || Wi) ) {
648: (*eps->ops->computevectors)(eps);
649: }
651: if (!eps->perm) k = i;
652: else k = eps->perm[i];
653: #ifdef PETSC_USE_COMPLEX
654: if (Wr) { VecCopy(eps->W[k], Wr); }
655: if (Wi) { VecSet(Wi,0.0); }
656: #else
657: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
658: if (Wr) { VecCopy(eps->W[k], Wr); }
659: if (Wi) { VecCopy(eps->W[k+1], Wi); }
660: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
661: if (Wr) { VecCopy(eps->W[k-1], Wr); }
662: if (Wi) {
663: VecCopy(eps->W[k], Wi);
664: VecScale(Wi,-1.0);
665: }
666: } else { /* real eigenvalue */
667: if (Wr) { VecCopy(eps->W[k], Wr); }
668: if (Wi) { VecSet(Wi,0.0); }
669: }
670: #endif
671:
672: return(0);
673: }
677: /*@
678: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
679: computed eigenpair.
681: Not Collective
683: Input Parameter:
684: + eps - eigensolver context
685: - i - index of eigenpair
687: Output Parameter:
688: . errest - the error estimate
690: Notes:
691: This is the error estimate used internally by the eigensolver. The actual
692: error bound can be computed with EPSComputeRelativeError(). See also the user's
693: manual for details.
695: Level: advanced
697: .seealso: EPSComputeRelativeError()
698: @*/
699: PetscErrorCode EPSGetErrorEstimate(EPS eps, PetscInt i, PetscReal *errest)
700: {
703: if (!eps->eigr || !eps->eigi) {
704: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
705: }
706: if (i<0 || i>=eps->nconv) {
707: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
708: }
709: if (eps->perm) i = eps->perm[i];
710: if (errest) *errest = eps->errest[i];
711: return(0);
712: }
716: /*@
717: EPSGetErrorEstimateLeft - Returns the left error estimate associated to the i-th
718: computed eigenpair (only available in two-sided eigensolvers).
720: Not Collective
722: Input Parameter:
723: + eps - eigensolver context
724: - i - index of eigenpair
726: Output Parameter:
727: . errest - the left error estimate
729: Notes:
730: This is the error estimate used internally by the eigensolver. The actual
731: error bound can be computed with EPSComputeRelativeErrorLeft(). See also the user's
732: manual for details.
734: Level: advanced
736: .seealso: EPSComputeRelativeErrorLeft()
737: @*/
738: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps, PetscInt i, PetscReal *errest)
739: {
742: if (!eps->eigr || !eps->eigi) {
743: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
744: }
745: if (eps->solverclass!=EPS_TWO_SIDE) {
746: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
747: }
748: if (i<0 || i>=eps->nconv) {
749: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
750: }
751: if (eps->perm) i = eps->perm[i];
752: if (errest) *errest = eps->errest_left[i];
753: return(0);
754: }
758: /*@
759: EPSComputeResidualNorm - Computes the norm of the residual vector associated with
760: the i-th computed eigenpair.
762: Collective on EPS
764: Input Parameter:
765: . eps - the eigensolver context
766: . i - the solution index
768: Output Parameter:
769: . norm - the residual norm, computed as ||Ax-kBx||_2 where k is the
770: eigenvalue and x is the eigenvector.
771: If k=0 then the residual norm is computed as ||Ax||_2.
773: Notes:
774: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
775: Eigenpairs are indexed according to the ordering criterion established
776: with EPSSetWhichEigenpairs().
778: Level: beginner
780: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
781: @*/
782: PetscErrorCode EPSComputeResidualNorm(EPS eps, PetscInt i, PetscReal *norm)
783: {
785: Vec u, v, w, xr, xi;
786: Mat A, B;
787: PetscScalar kr, ki;
788: #ifndef PETSC_USE_COMPLEX
789: PetscReal ni, nr;
790: #endif
791:
794: STGetOperators(eps->OP,&A,&B);
795: VecDuplicate(eps->vec_initial,&u);
796: VecDuplicate(eps->vec_initial,&v);
797: VecDuplicate(eps->vec_initial,&w);
798: VecDuplicate(eps->vec_initial,&xr);
799: VecDuplicate(eps->vec_initial,&xi);
800: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
802: #ifndef PETSC_USE_COMPLEX
803: if (ki == 0 ||
804: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
805: #endif
806: MatMult( A, xr, u ); /* u=A*x */
807: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
808: if (eps->isgeneralized) { MatMult( B, xr, w ); }
809: else { VecCopy( xr, w ); } /* w=B*x */
810: VecAXPY( u, -kr, w ); /* u=A*x-k*B*x */
811: }
812: VecNorm( u, NORM_2, norm);
813: #ifndef PETSC_USE_COMPLEX
814: } else {
815: MatMult( A, xr, u ); /* u=A*xr */
816: if (eps->isgeneralized) { MatMult( B, xr, v ); }
817: else { VecCopy( xr, v ); } /* v=B*xr */
818: VecAXPY( u, -kr, v ); /* u=A*xr-kr*B*xr */
819: if (eps->isgeneralized) { MatMult( B, xi, w ); }
820: else { VecCopy( xi, w ); } /* w=B*xi */
821: VecAXPY( u, ki, w ); /* u=A*xr-kr*B*xr+ki*B*xi */
822: VecNorm( u, NORM_2, &nr );
823: MatMult( A, xi, u ); /* u=A*xi */
824: VecAXPY( u, -kr, w ); /* u=A*xi-kr*B*xi */
825: VecAXPY( u, -ki, v ); /* u=A*xi-kr*B*xi-ki*B*xr */
826: VecNorm( u, NORM_2, &ni );
827: *norm = SlepcAbsEigenvalue( nr, ni );
828: }
829: #endif
831: VecDestroy(w);
832: VecDestroy(v);
833: VecDestroy(u);
834: VecDestroy(xr);
835: VecDestroy(xi);
836: return(0);
837: }
841: /*@
842: EPSComputeResidualNormLeft - Computes the norm of the residual vector associated with
843: the i-th computed left eigenvector (only available in two-sided eigensolvers).
845: Collective on EPS
847: Input Parameter:
848: . eps - the eigensolver context
849: . i - the solution index
851: Output Parameter:
852: . norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the
853: eigenvalue and y is the left eigenvector.
854: If k=0 then the residual norm is computed as ||y'A||_2.
856: Notes:
857: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
858: Eigenpairs are indexed according to the ordering criterion established
859: with EPSSetWhichEigenpairs().
861: Level: beginner
863: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
864: @*/
865: PetscErrorCode EPSComputeResidualNormLeft(EPS eps, PetscInt i, PetscReal *norm)
866: {
868: Vec u, v, w, xr, xi;
869: Mat A, B;
870: PetscScalar kr, ki;
871: #ifndef PETSC_USE_COMPLEX
872: PetscReal ni, nr;
873: #endif
874:
877: STGetOperators(eps->OP,&A,&B);
878: VecDuplicate(eps->vec_initial_left,&u);
879: VecDuplicate(eps->vec_initial_left,&v);
880: VecDuplicate(eps->vec_initial_left,&w);
881: VecDuplicate(eps->vec_initial_left,&xr);
882: VecDuplicate(eps->vec_initial_left,&xi);
883: EPSGetValue(eps,i,&kr,&ki);
884: EPSGetLeftVector(eps,i,xr,xi);
886: #ifndef PETSC_USE_COMPLEX
887: if (ki == 0 ||
888: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
889: #endif
890: MatMultTranspose( A, xr, u ); /* u=A'*x */
891: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
892: if (eps->isgeneralized) { MatMultTranspose( B, xr, w ); }
893: else { VecCopy( xr, w ); } /* w=B'*x */
894: VecAXPY( u, -kr, w); /* u=A'*x-k*B'*x */
895: }
896: VecNorm( u, NORM_2, norm);
897: #ifndef PETSC_USE_COMPLEX
898: } else {
899: MatMultTranspose( A, xr, u ); /* u=A'*xr */
900: if (eps->isgeneralized) { MatMultTranspose( B, xr, v ); }
901: else { VecCopy( xr, v ); } /* v=B'*xr */
902: VecAXPY( u, -kr, v ); /* u=A'*xr-kr*B'*xr */
903: if (eps->isgeneralized) { MatMultTranspose( B, xi, w ); }
904: else { VecCopy( xi, w ); } /* w=B'*xi */
905: VecAXPY( u, ki, w ); /* u=A'*xr-kr*B'*xr+ki*B'*xi */
906: VecNorm( u, NORM_2, &nr );
907: MatMultTranspose( A, xi, u ); /* u=A'*xi */
908: VecAXPY( u, -kr, w ); /* u=A'*xi-kr*B'*xi */
909: VecAXPY( u, -ki, v ); /* u=A'*xi-kr*B'*xi-ki*B'*xr */
910: VecNorm( u, NORM_2, &ni );
911: *norm = SlepcAbsEigenvalue( nr, ni );
912: }
913: #endif
915: VecDestroy(w);
916: VecDestroy(v);
917: VecDestroy(u);
918: VecDestroy(xr);
919: VecDestroy(xi);
920: return(0);
921: }
925: /*@
926: EPSComputeRelativeError - Computes the relative error bound associated
927: with the i-th computed eigenpair.
929: Collective on EPS
931: Input Parameter:
932: . eps - the eigensolver context
933: . i - the solution index
935: Output Parameter:
936: . error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where
937: k is the eigenvalue and x is the eigenvector.
938: If k=0 the relative error is computed as ||Ax||_2/||x||_2.
940: Level: beginner
942: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
943: @*/
944: PetscErrorCode EPSComputeRelativeError(EPS eps, PetscInt i, PetscReal *error)
945: {
947: Vec xr, xi;
948: PetscScalar kr, ki;
949: PetscReal norm, er;
950: #ifndef PETSC_USE_COMPLEX
951: Vec u;
952: PetscReal ei;
953: #endif
954:
957: EPSComputeResidualNorm(eps,i,&norm);
958: VecDuplicate(eps->vec_initial,&xr);
959: VecDuplicate(eps->vec_initial,&xi);
960: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
962: #ifndef PETSC_USE_COMPLEX
963: if (ki == 0 ||
964: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
965: #endif
966: VecNorm(xr, NORM_2, &er);
967: if (PetscAbsScalar(kr) > norm) {
968: *error = norm / (PetscAbsScalar(kr) * er);
969: } else {
970: *error = norm / er;
971: }
972: #ifndef PETSC_USE_COMPLEX
973: } else {
974: if (SlepcAbsEigenvalue(kr,ki) > norm) {
975: VecDuplicate(xi, &u);
976: VecCopy(xi, u);
977: VecAXPBY(u, kr, -ki, xr);
978: VecNorm(u, NORM_2, &er);
979: VecAXPBY(xi, kr, ki, xr);
980: VecNorm(xi, NORM_2, &ei);
981: VecDestroy(u);
982: } else {
983: VecDot(xr, xr, &er);
984: VecDot(xi, xi, &ei);
985: }
986: *error = norm / SlepcAbsEigenvalue(er, ei);
987: }
988: #endif
989:
990: VecDestroy(xr);
991: VecDestroy(xi);
992: return(0);
993: }
997: /*@
998: EPSComputeRelativeErrorLeft - Computes the relative error bound associated
999: with the i-th computed eigenvalue and left eigenvector (only available in
1000: two-sided eigensolvers).
1002: Collective on EPS
1004: Input Parameter:
1005: . eps - the eigensolver context
1006: . i - the solution index
1008: Output Parameter:
1009: . error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where
1010: k is the eigenvalue and y is the left eigenvector.
1011: If k=0 the relative error is computed as ||y'A||_2/||y||_2.
1013: Level: beginner
1015: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
1016: @*/
1017: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps, PetscInt i, PetscReal *error)
1018: {
1020: Vec xr, xi;
1021: PetscScalar kr, ki;
1022: PetscReal norm, er;
1023: #ifndef PETSC_USE_COMPLEX
1024: Vec u;
1025: PetscReal ei;
1026: #endif
1027:
1030: EPSComputeResidualNormLeft(eps,i,&norm);
1031: VecDuplicate(eps->vec_initial_left,&xr);
1032: VecDuplicate(eps->vec_initial_left,&xi);
1033: EPSGetValue(eps,i,&kr,&ki);
1034: EPSGetLeftVector(eps,i,xr,xi);
1036: #ifndef PETSC_USE_COMPLEX
1037: if (ki == 0 ||
1038: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1039: #endif
1040: VecNorm(xr, NORM_2, &er);
1041: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1042: *error = norm / (PetscAbsScalar(kr) * er);
1043: } else {
1044: *error = norm / er;
1045: }
1046: #ifndef PETSC_USE_COMPLEX
1047: } else {
1048: VecDuplicate(xi, &u);
1049: VecCopy(xi, u);
1050: VecAXPBY(u, kr, -ki, xr);
1051: VecNorm(u, NORM_2, &er);
1052: VecAXPBY(xi, kr, ki, xr);
1053: VecNorm(xi, NORM_2, &ei);
1054: VecDestroy(u);
1055: *error = norm / SlepcAbsEigenvalue(er, ei);
1056: }
1057: #endif
1058:
1059: VecDestroy(xr);
1060: VecDestroy(xi);
1061: return(0);
1062: }
1064: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1068: /*@
1069: EPSSortEigenvalues - Sorts a list of eigenvalues according to a certain
1070: criterion.
1072: Not Collective
1074: Input Parameters:
1075: + n - number of eigenvalue in the list
1076: . eig - pointer to the array containing the eigenvalues
1077: . eigi - imaginary part of the eigenvalues (only when using real numbers)
1078: . which - sorting criterion
1079: - nev - number of wanted eigenvalues
1081: Output Parameter:
1082: . permout - resulting permutation
1084: Notes:
1085: The result is a list of indices in the original eigenvalue array
1086: corresponding to the first nev eigenvalues sorted in the specified
1087: criterion
1089: Level: developer
1091: .seealso: EPSSortEigenvaluesReal(), EPSSetWhichEigenpairs()
1092: @*/
1093: PetscErrorCode EPSSortEigenvalues(PetscInt n,PetscScalar *eig,PetscScalar *eigi,EPSWhich which,PetscInt nev,PetscInt *permout)
1094: {
1096: PetscInt i;
1097: PetscInt *perm;
1098: PetscReal *values;
1101: PetscMalloc(n*sizeof(PetscInt),&perm);
1102: PetscMalloc(n*sizeof(PetscReal),&values);
1103: for (i=0; i<n; i++) { perm[i] = i;}
1105: switch(which) {
1106: case EPS_LARGEST_MAGNITUDE:
1107: case EPS_SMALLEST_MAGNITUDE:
1108: for (i=0; i<n; i++) { values[i] = SlepcAbsEigenvalue(eig[i],eigi[i]); }
1109: break;
1110: case EPS_LARGEST_REAL:
1111: case EPS_SMALLEST_REAL:
1112: for (i=0; i<n; i++) { values[i] = PetscRealPart(eig[i]); }
1113: break;
1114: case EPS_LARGEST_IMAGINARY:
1115: case EPS_SMALLEST_IMAGINARY:
1116: #if defined(PETSC_USE_COMPLEX)
1117: for (i=0; i<n; i++) { values[i] = PetscImaginaryPart(eig[i]); }
1118: #else
1119: for (i=0; i<n; i++) { values[i] = PetscAbsReal(eigi[i]); }
1120: #endif
1121: break;
1122: default: SETERRQ(1,"Wrong value of which");
1123: }
1125: PetscSortRealWithPermutation(n,values,perm);
1127: switch(which) {
1128: case EPS_LARGEST_MAGNITUDE:
1129: case EPS_LARGEST_REAL:
1130: case EPS_LARGEST_IMAGINARY:
1131: for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
1132: break;
1133: case EPS_SMALLEST_MAGNITUDE:
1134: case EPS_SMALLEST_REAL:
1135: case EPS_SMALLEST_IMAGINARY:
1136: for (i=0; i<nev; i++) { permout[i] = perm[i]; }
1137: break;
1138: default: SETERRQ(1,"Wrong value of which");
1139: }
1141: #if !defined(PETSC_USE_COMPLEX)
1142: for (i=0; i<nev-1; i++) {
1143: if (eigi[permout[i]] != 0.0) {
1144: if (eig[permout[i]] == eig[permout[i+1]] &&
1145: eigi[permout[i]] == -eigi[permout[i+1]] &&
1146: eigi[permout[i]] < 0.0) {
1147: PetscInt tmp;
1148: SWAP(permout[i], permout[i+1], tmp);
1149: }
1150: i++;
1151: }
1152: }
1153: #endif
1155: PetscFree(values);
1156: PetscFree(perm);
1157: return(0);
1158: }
1162: /*@
1163: EPSSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
1164: criterion (version for real eigenvalues only).
1166: Not Collective
1168: Input Parameters:
1169: + n - number of eigenvalue in the list
1170: . eig - pointer to the array containing the eigenvalues (real)
1171: . which - sorting criterion
1172: - nev - number of wanted eigenvalues
1174: Output Parameter:
1175: . permout - resulting permutation
1177: Workspace:
1178: . work - workspace for storing n real values and n integer values
1180: Notes:
1181: The result is a list of indices in the original eigenvalue array
1182: corresponding to the first nev eigenvalues sorted in the specified
1183: criterion
1185: Level: developer
1187: .seealso: EPSSortEigenvalues(), EPSSetWhichEigenpairs()
1188: @*/
1189: PetscErrorCode EPSSortEigenvaluesReal(PetscInt n,PetscReal *eig,EPSWhich which,PetscInt nev,PetscInt *permout,PetscReal *work)
1190: {
1192: PetscInt i;
1193: PetscReal *values = work;
1194: PetscInt *perm = (PetscInt*)(work+n);
1197: for (i=0; i<n; i++) { perm[i] = i;}
1199: switch(which) {
1200: case EPS_LARGEST_MAGNITUDE:
1201: case EPS_SMALLEST_MAGNITUDE:
1202: for (i=0; i<n; i++) { values[i] = PetscAbsReal(eig[i]); }
1203: break;
1204: case EPS_LARGEST_REAL:
1205: case EPS_SMALLEST_REAL:
1206: for (i=0; i<n; i++) { values[i] = eig[i]; }
1207: break;
1208: default: SETERRQ(1,"Wrong value of which");
1209: }
1211: PetscSortRealWithPermutation(n,values,perm);
1213: switch(which) {
1214: case EPS_LARGEST_MAGNITUDE:
1215: case EPS_LARGEST_REAL:
1216: for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
1217: break;
1218: case EPS_SMALLEST_MAGNITUDE:
1219: case EPS_SMALLEST_REAL:
1220: for (i=0; i<nev; i++) { permout[i] = perm[i]; }
1221: break;
1222: default: SETERRQ(1,"Wrong value of which");
1223: }
1224: return(0);
1225: }
1229: /*@
1230: EPSGetStartVector - Gets a vector to be used as the starting vector
1231: in an Arnoldi or Lanczos reduction.
1233: Collective on EPS and Vec
1235: Input Parameters:
1236: + eps - the eigensolver context
1237: - i - index of the Arnoldi/Lanczos step
1239: Output Parameters:
1240: + vec - the start vector
1241: - breakdown - flag indicating that a breakdown has occurred
1243: Notes:
1244: The start vector is computed from another vector: for the first step (i=0),
1245: the initial vector is used (see EPSGetInitialVector()); otherwise a random
1246: vector is created. Then this vector is forced to be in the range of OP (only
1247: for generalized definite problems) and orthonormalized with respect to all
1248: V-vectors up to i-1.
1250: The flag breakdown is set to true if either i=0 and the vector belongs to the
1251: deflation space, or i>0 and the vector is linearly dependent with respect
1252: to the V-vectors.
1254: The caller must pass a vector already allocated with dimensions conforming
1255: to the initial vector. This vector is overwritten.
1257: Level: developer
1259: .seealso: EPSGetInitialVector()
1261: @*/
1262: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,Vec vec,PetscTruth *breakdown)
1263: {
1265: PetscReal norm;
1266: PetscTruth lindep;
1267: Vec w;
1268:
1273: /* For the first step, use the initial vector, otherwise a random one */
1274: if (i==0) {
1275: w = eps->vec_initial;
1276: } else {
1277: VecDuplicate(eps->vec_initial,&w);
1278: SlepcVecSetRandom(w);
1279: }
1281: /* Force the vector to be in the range of OP for definite generalized problems */
1282: if (eps->ispositive) {
1283: STApply(eps->OP,w,vec);
1284: } else {
1285: VecCopy(w,vec);
1286: }
1288: /* Orthonormalize the vector with respect to previous vectors */
1289: IPOrthogonalize(eps->ip,i+eps->nds,PETSC_NULL,eps->DSV,vec,PETSC_NULL,&norm,&lindep,PETSC_NULL,PETSC_NULL);
1290: if (breakdown) *breakdown = lindep;
1291: else if (lindep || norm == 0.0) {
1292: if (i==0) { SETERRQ(1,"Initial vector is zero or belongs to the deflation space"); }
1293: else { SETERRQ(1,"Unable to generate more start vectors"); }
1294: }
1295:
1296: VecScale(vec,1.0/norm);
1298: if (i!=0) {
1299: VecDestroy(w);
1300: }
1302: return(0);
1303: }
1306: /*@
1307: EPSGetLeftStartVector - Gets a vector to be used as the starting vector
1308: in the left recurrence of a two-sided eigensolver.
1310: Collective on EPS and Vec
1312: Input Parameters:
1313: + eps - the eigensolver context
1314: - i - index of the Arnoldi/Lanczos step
1316: Output Parameter:
1317: . vec - the start vector
1319: Notes:
1320: The start vector is computed from another vector: for the first step (i=0),
1321: the left initial vector is used (see EPSGetLeftInitialVector()); otherwise
1322: a random vector is created. Then this vector is forced to be in the range
1323: of OP' and orthonormalized with respect to all W-vectors up to i-1.
1325: The caller must pass a vector already allocated with dimensions conforming
1326: to the left initial vector. This vector is overwritten.
1328: Level: developer
1330: .seealso: EPSGetLeftInitialVector()
1332: @*/
1333: PetscErrorCode EPSGetLeftStartVector(EPS eps,PetscInt i,Vec vec)
1334: {
1336: PetscTruth breakdown;
1337: PetscReal norm;
1338: Vec w;
1339:
1344: /* For the first step, use the initial vector, otherwise a random one */
1345: if (i==0) {
1346: w = eps->vec_initial_left;
1347: }
1348: else {
1349: VecDuplicate(eps->vec_initial_left,&w);
1350: SlepcVecSetRandom(w);
1351: }
1353: /* Force the vector to be in the range of OP */
1354: STApplyTranspose(eps->OP,w,vec);
1356: /* Orthonormalize the vector with respect to previous vectors */
1357: IPOrthogonalize(eps->ip,i,PETSC_NULL,eps->W,vec,PETSC_NULL,&norm,&breakdown,PETSC_NULL,PETSC_NULL);
1358: if (breakdown) {
1359: if (i==0) { SETERRQ(1,"Left initial vector is zero"); }
1360: else { SETERRQ(1,"Unable to generate more left start vectors"); }
1361: }
1362: VecScale(vec,1/norm);
1364: if (i!=0) {
1365: VecDestroy(w);
1366: }
1368: return(0);
1369: }