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-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/eps/epsimpl.h
17: /*@
18: EPSSolve - Solves the eigensystem.
20: Collective on EPS
22: Input Parameter:
23: . eps - eigensolver context obtained from EPSCreate()
25: Options Database:
26: + -eps_view - print information about the solver used
27: . -eps_view_binary - save the matrices to the default binary file
28: - -eps_plot_eigs - plot computed eigenvalues
30: Level: beginner
32: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
33: @*/
34: PetscErrorCode EPSSolve(EPS eps)
35: {
37: int i;
38: PetscReal re,im;
39: PetscTruth flg;
40: PetscViewer viewer;
41: PetscDraw draw;
42: PetscDrawSP drawsp;
47: PetscOptionsHasName(eps->prefix,"-eps_view_binary",&flg);
48: if (flg) {
49: Mat A,B;
50: STGetOperators(eps->OP,&A,&B);
51: MatView(A,PETSC_VIEWER_BINARY_(eps->comm));
52: if (B) MatView(B,PETSC_VIEWER_BINARY_(eps->comm));
53: }
55: /* reset the convergence flag from the previous solves */
56: eps->reason = EPS_CONVERGED_ITERATING;
58: if (!eps->setupcalled){ EPSSetUp(eps); }
59: STResetOperationCounters(eps->OP);
60: eps->nv = eps->ncv;
61: eps->evecsavailable = PETSC_FALSE;
62: eps->nconv = 0;
63: eps->its = 0;
64: for (i=0;i<eps->ncv;i++) eps->eigr[i]=eps->eigi[i]=eps->errest[i]=0.0;
65: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nv);
67: PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
69: switch (eps->solverclass) {
70: case EPS_ONE_SIDE:
71: (*eps->ops->solve)(eps); break;
72: case EPS_TWO_SIDE:
73: if (eps->ops->solvets) {
74: (*eps->ops->solvets)(eps); break;
75: } else {
76: SETERRQ(1,"Two-sided version unavailable for this solver");
77: }
78: default:
79: SETERRQ(1,"Wrong value of eps->solverclass");
80: }
82: STPostSolve(eps->OP);
83: PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
85: if (!eps->reason) {
86: SETERRQ(1,"Internal error, solver returned without setting converged reason");
87: }
89: /* Map eigenvalues back to the original problem, necessary in some
90: * spectral transformations */
91: (*eps->ops->backtransform)(eps);
93: /* Adjust left eigenvectors in generalized problems: y = B^T y */
94: if (eps->isgeneralized && eps->solverclass == EPS_TWO_SIDE) {
95: Mat B;
96: KSP ksp;
97: Vec w;
98: STGetOperators(eps->OP,PETSC_NULL,&B);
99: KSPCreate(eps->comm,&ksp);
100: KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
101: KSPSetFromOptions(ksp);
102: MatGetVecs(B,PETSC_NULL,&w);
103: for (i=0;i<eps->nconv;i++) {
104: VecCopy(eps->W[i],w);
105: KSPSolveTranspose(ksp,w,eps->W[i]);
106: }
107: KSPDestroy(ksp);
108: VecDestroy(w);
109: }
110:
112: #ifndef PETSC_USE_COMPLEX
113: /* reorder conjugate eigenvalues (positive imaginary first) */
114: for (i=0; i<eps->nconv-1; i++) {
115: if (eps->eigi[i] != 0) {
116: if (eps->eigi[i] < 0) {
117: eps->eigi[i] = -eps->eigi[i];
118: eps->eigi[i+1] = -eps->eigi[i+1];
119: VecScale(eps->V[i+1],-1.0);
120: }
121: i++;
122: }
123: }
124: #endif
126: /* sort eigenvalues according to eps->which parameter */
127: PetscFree(eps->perm);
128: if (eps->nconv > 0) {
129: PetscMalloc(sizeof(int)*eps->nconv, &eps->perm);
130: EPSSortEigenvalues(eps->nconv, eps->eigr, eps->eigi, eps->which, eps->nconv, eps->perm);
131: }
133: PetscOptionsHasName(eps->prefix,"-eps_view",&flg);
134: if (flg && !PetscPreLoadingOn) { EPSView(eps,PETSC_VIEWER_STDOUT_WORLD); }
136: PetscOptionsHasName(eps->prefix,"-eps_plot_eigs",&flg);
137: if (flg) {
138: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
139: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
140: PetscViewerDrawGetDraw(viewer,0,&draw);
141: PetscDrawSPCreate(draw,1,&drawsp);
142: for( i=0; i<eps->nconv; i++ ) {
143: #if defined(PETSC_USE_COMPLEX)
144: re = PetscRealPart(eps->eigr[i]);
145: im = PetscImaginaryPart(eps->eigi[i]);
146: #else
147: re = eps->eigr[i];
148: im = eps->eigi[i];
149: #endif
150: PetscDrawSPAddPoint(drawsp,&re,&im);
151: }
152: PetscDrawSPDraw(drawsp);
153: PetscDrawSPDestroy(drawsp);
154: PetscViewerDestroy(viewer);
155: }
157: return(0);
158: }
162: /*@
163: EPSGetIterationNumber - Gets the current iteration number. If the
164: call to EPSSolve() is complete, then it returns the number of iterations
165: carried out by the solution method.
166:
167: Not Collective
169: Input Parameter:
170: . eps - the eigensolver context
172: Output Parameter:
173: . its - number of iterations
175: Level: intermediate
177: Note:
178: During the i-th iteration this call returns i-1. If EPSSolve() is
179: complete, then parameter "its" contains either the iteration number at
180: which convergence was successfully reached, or failure was detected.
181: Call EPSGetConvergedReason() to determine if the solver converged or
182: failed and why.
184: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
185: @*/
186: PetscErrorCode EPSGetIterationNumber(EPS eps,int *its)
187: {
191: *its = eps->its;
192: return(0);
193: }
197: /*@
198: EPSGetOperationCounters - Gets the total number of operator applications,
199: inner product operations and linear iterations used by the ST object
200: during the last EPSSolve() call.
202: Not Collective
204: Input Parameter:
205: . eps - EPS context
207: Output Parameter:
208: + ops - number of operator applications
209: . dots - number of inner product operations
210: - lits - number of linear iterations
212: Notes:
213: When the eigensolver algorithm invokes STApply() then a linear system
214: must be solved (except in the case of standard eigenproblems and shift
215: transformation). The number of iterations required in this solve is
216: accumulated into a counter whose value is returned by this function.
218: These counters are reset to zero at each successive call to EPSSolve().
220: Level: intermediate
222: @*/
223: PetscErrorCode EPSGetOperationCounters(EPS eps,int* ops,int* dots,int* lits)
224: {
229: STGetOperationCounters(eps->OP,ops,lits);
230: if (dots) {
231: IPGetOperationCounters(eps->ip,dots);
232: }
233: return(0);
234: }
238: /*@
239: EPSGetConverged - Gets the number of converged eigenpairs.
241: Not Collective
243: Input Parameter:
244: . eps - the eigensolver context
245:
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()
255: @*/
256: PetscErrorCode EPSGetConverged(EPS eps,int *nconv)
257: {
261: *nconv = eps->nconv;
262: return(0);
263: }
268: /*@C
269: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
270: stopped.
272: Not Collective
274: Input Parameter:
275: . eps - the eigensolver context
277: Output Parameter:
278: . reason - negative value indicates diverged, positive value converged
279: (see EPSConvergedReason)
281: Possible values for reason:
282: + EPS_CONVERGED_TOL - converged up to tolerance
283: . EPS_DIVERGED_ITS - required more than its to reach convergence
284: . EPS_DIVERGED_BREAKDOWN - generic breakdown in method
285: - EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric
287: Level: intermediate
289: Notes: Can only be called after the call to EPSSolve() is complete.
291: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
292: @*/
293: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
294: {
298: *reason = eps->reason;
299: return(0);
300: }
304: /*@
305: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
306: subspace.
308: Not Collective
310: Input Parameter:
311: . eps - the eigensolver context
312:
313: Output Parameter:
314: . v - an array of vectors
316: Notes:
317: This function should be called after EPSSolve() has finished.
319: The user should provide in v an array of nconv vectors, where nconv is
320: the value returned by EPSGetConverged().
322: The first k vectors returned in v span an invariant subspace associated
323: with the first k computed eigenvalues (note that this is not true if the
324: k-th eigenvalue is complex and matrix A is real; in this case the first
325: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
326: in X for all x in X (a similar definition applies for generalized
327: eigenproblems).
329: Level: intermediate
331: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetLeftInvariantSubspace()
332: @*/
333: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
334: {
336: int i;
342: if (!eps->V) {
343: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
344: }
345: for (i=0;i<eps->nconv;i++) {
346: VecCopy(eps->V[i],v[i]);
347: }
348: return(0);
349: }
353: /*@
354: EPSGetLeftInvariantSubspace - Gets an orthonormal basis of the computed left
355: invariant subspace (only available in two-sided eigensolvers).
357: Not Collective
359: Input Parameter:
360: . eps - the eigensolver context
361:
362: Output Parameter:
363: . v - an array of vectors
365: Notes:
366: This function should be called after EPSSolve() has finished.
368: The user should provide in v an array of nconv vectors, where nconv is
369: the value returned by EPSGetConverged().
371: The first k vectors returned in v span a left invariant subspace associated
372: with the first k computed eigenvalues (note that this is not true if the
373: k-th eigenvalue is complex and matrix A is real; in this case the first
374: k+1 vectors should be used). A left invariant subspace Y of A satisfies y'A
375: in Y for all y in Y (a similar definition applies for generalized
376: eigenproblems).
378: Level: intermediate
380: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
381: @*/
382: PetscErrorCode EPSGetLeftInvariantSubspace(EPS eps, Vec *v)
383: {
385: int i;
391: if (!eps->W) {
392: if (eps->solverclass!=EPS_TWO_SIDE) {
393: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
394: } else {
395: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
396: }
397: }
398: for (i=0;i<eps->nconv;i++) {
399: VecCopy(eps->W[i],v[i]);
400: }
401: return(0);
402: }
406: /*@
407: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
408: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
410: Not Collective
412: Input Parameters:
413: + eps - eigensolver context
414: - i - index of the solution
416: Output Parameters:
417: + eigr - real part of eigenvalue
418: . eigi - imaginary part of eigenvalue
419: . Vr - real part of eigenvector
420: - Vi - imaginary part of eigenvector
422: Notes:
423: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
424: configured with complex scalars the eigenvalue is stored
425: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
426: set to zero).
428: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
429: Eigenpairs are indexed according to the ordering criterion established
430: with EPSSetWhichEigenpairs().
432: Level: beginner
434: .seealso: EPSGetValue(), EPSGetRightVector(), EPSGetLeftVector(), EPSSolve(),
435: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
436: @*/
437: PetscErrorCode EPSGetEigenpair(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
438: {
443: if (!eps->eigr || !eps->eigi || !eps->V) {
444: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
445: }
446: if (i<0 || i>=eps->nconv) {
447: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
448: }
449: EPSGetValue(eps,i,eigr,eigi);
450: EPSGetRightVector(eps,i,Vr,Vi);
451:
452: return(0);
453: }
457: /*@
458: EPSGetValue - Gets the i-th eigenvalue as computed by EPSSolve().
460: Not Collective
462: Input Parameters:
463: + eps - eigensolver context
464: - i - index of the solution
466: Output Parameters:
467: + eigr - real part of eigenvalue
468: - eigi - imaginary part of eigenvalue
470: Notes:
471: If the eigenvalue is real, then eigi is set to zero. If PETSc is
472: configured with complex scalars the eigenvalue is stored
473: directly in eigr (eigi is set to zero).
475: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
476: Eigenpairs are indexed according to the ordering criterion established
477: with EPSSetWhichEigenpairs().
479: Level: beginner
481: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
482: EPSGetEigenpair()
483: @*/
484: PetscErrorCode EPSGetValue(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi)
485: {
486: int k;
490: if (!eps->eigr || !eps->eigi) {
491: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
492: }
493: if (i<0 || i>=eps->nconv) {
494: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
495: }
497: if (!eps->perm) k = i;
498: else k = eps->perm[i];
499: #ifdef PETSC_USE_COMPLEX
500: if (eigr) *eigr = eps->eigr[k];
501: if (eigi) *eigi = 0;
502: #else
503: if (eigr) *eigr = eps->eigr[k];
504: if (eigi) *eigi = eps->eigi[k];
505: #endif
506:
507: return(0);
508: }
512: /*@
513: EPSGetRightVector - Gets the i-th right eigenvector as computed by EPSSolve().
515: Not Collective
517: Input Parameters:
518: + eps - eigensolver context
519: - i - index of the solution
521: Output Parameters:
522: + Vr - real part of eigenvector
523: - Vi - imaginary part of eigenvector
525: Notes:
526: If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
527: configured with complex scalars the eigenvector is stored
528: directly in Vr (Vi is set to zero).
530: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
531: Eigenpairs are indexed according to the ordering criterion established
532: with EPSSetWhichEigenpairs().
534: Level: beginner
536: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
537: EPSGetEigenpair(), EPSGetLeftVector()
538: @*/
539: PetscErrorCode EPSGetRightVector(EPS eps, int i, Vec Vr, Vec Vi)
540: {
542: int k;
546: if (!eps->V) {
547: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
548: }
549: if (i<0 || i>=eps->nconv) {
550: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
551: }
552: if (!eps->evecsavailable && (Vr || Vi) ) {
553: (*eps->ops->computevectors)(eps);
554: }
556: if (!eps->perm) k = i;
557: else k = eps->perm[i];
558: #ifdef PETSC_USE_COMPLEX
559: if (Vr) { VecCopy(eps->AV[k], Vr); }
560: if (Vi) { VecSet(Vi,0.0); }
561: #else
562: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
563: if (Vr) { VecCopy(eps->AV[k], Vr); }
564: if (Vi) { VecCopy(eps->AV[k+1], Vi); }
565: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
566: if (Vr) { VecCopy(eps->AV[k-1], Vr); }
567: if (Vi) {
568: VecCopy(eps->AV[k], Vi);
569: VecScale(Vi,-1.0);
570: }
571: } else { /* real eigenvalue */
572: if (Vr) { VecCopy(eps->AV[k], Vr); }
573: if (Vi) { VecSet(Vi,0.0); }
574: }
575: #endif
576:
577: return(0);
578: }
582: /*@
583: EPSGetLeftVector - Gets the i-th left eigenvector as computed by EPSSolve()
584: (only available in two-sided eigensolvers).
586: Not Collective
588: Input Parameters:
589: + eps - eigensolver context
590: - i - index of the solution
592: Output Parameters:
593: + Wr - real part of eigenvector
594: - Wi - imaginary part of eigenvector
596: Notes:
597: If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is
598: configured with complex scalars the eigenvector is stored
599: directly in Wr (Wi is set to zero).
601: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
602: Eigenpairs are indexed according to the ordering criterion established
603: with EPSSetWhichEigenpairs().
605: Level: beginner
607: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
608: EPSGetEigenpair(), EPSGetLeftVector()
609: @*/
610: PetscErrorCode EPSGetLeftVector(EPS eps, int i, Vec Wr, Vec Wi)
611: {
613: int k;
617: if (!eps->W) {
618: if (eps->solverclass!=EPS_TWO_SIDE) {
619: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
620: } else {
621: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
622: }
623: }
624: if (i<0 || i>=eps->nconv) {
625: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
626: }
627: if (!eps->evecsavailable && (Wr || Wi) ) {
628: (*eps->ops->computevectors)(eps);
629: }
631: if (!eps->perm) k = i;
632: else k = eps->perm[i];
633: #ifdef PETSC_USE_COMPLEX
634: if (Wr) { VecCopy(eps->AW[k], Wr); }
635: if (Wi) { VecSet(Wi,0.0); }
636: #else
637: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
638: if (Wr) { VecCopy(eps->AW[k], Wr); }
639: if (Wi) { VecCopy(eps->AW[k+1], Wi); }
640: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
641: if (Wr) { VecCopy(eps->AW[k-1], Wr); }
642: if (Wi) {
643: VecCopy(eps->AW[k], Wi);
644: VecScale(Wi,-1.0);
645: }
646: } else { /* real eigenvalue */
647: if (Wr) { VecCopy(eps->AW[k], Wr); }
648: if (Wi) { VecSet(Wi,0.0); }
649: }
650: #endif
651:
652: return(0);
653: }
657: /*@
658: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
659: computed eigenpair.
661: Not Collective
663: Input Parameter:
664: + eps - eigensolver context
665: - i - index of eigenpair
667: Output Parameter:
668: . errest - the error estimate
670: Notes:
671: This is the error estimate used internally by the eigensolver. The actual
672: error bound can be computed with EPSComputeRelativeError(). See also the user's
673: manual for details.
675: Level: advanced
677: .seealso: EPSComputeRelativeError()
678: @*/
679: PetscErrorCode EPSGetErrorEstimate(EPS eps, int i, PetscReal *errest)
680: {
683: if (!eps->eigr || !eps->eigi) {
684: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
685: }
686: if (i<0 || i>=eps->nconv) {
687: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
688: }
689: if (eps->perm) i = eps->perm[i];
690: if (errest) *errest = eps->errest[i];
691: return(0);
692: }
696: /*@
697: EPSGetErrorEstimateLeft - Returns the left error estimate associated to the i-th
698: computed eigenpair (only available in two-sided eigensolvers).
700: Not Collective
702: Input Parameter:
703: + eps - eigensolver context
704: - i - index of eigenpair
706: Output Parameter:
707: . errest - the left error estimate
709: Notes:
710: This is the error estimate used internally by the eigensolver. The actual
711: error bound can be computed with EPSComputeRelativeErrorLeft(). See also the user's
712: manual for details.
714: Level: advanced
716: .seealso: EPSComputeRelativeErrorLeft()
717: @*/
718: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps, int i, PetscReal *errest)
719: {
722: if (!eps->eigr || !eps->eigi) {
723: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
724: }
725: if (eps->solverclass!=EPS_TWO_SIDE) {
726: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
727: }
728: if (i<0 || i>=eps->nconv) {
729: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
730: }
731: if (eps->perm) i = eps->perm[i];
732: if (errest) *errest = eps->errest_left[i];
733: return(0);
734: }
738: /*@
739: EPSComputeResidualNorm - Computes the norm of the residual vector associated with
740: the i-th computed eigenpair.
742: Collective on EPS
744: Input Parameter:
745: . eps - the eigensolver context
746: . i - the solution index
748: Output Parameter:
749: . norm - the residual norm, computed as ||Ax-kBx||_2 where k is the
750: eigenvalue and x is the eigenvector.
751: If k=0 then the residual norm is computed as ||Ax||_2.
753: Notes:
754: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
755: Eigenpairs are indexed according to the ordering criterion established
756: with EPSSetWhichEigenpairs().
758: Level: beginner
760: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
761: @*/
762: PetscErrorCode EPSComputeResidualNorm(EPS eps, int i, PetscReal *norm)
763: {
765: Vec u, v, w, xr, xi;
766: Mat A, B;
767: PetscScalar kr, ki;
768: #ifndef PETSC_USE_COMPLEX
769: PetscReal ni, nr;
770: #endif
771:
774: STGetOperators(eps->OP,&A,&B);
775: VecDuplicate(eps->vec_initial,&u);
776: VecDuplicate(eps->vec_initial,&v);
777: VecDuplicate(eps->vec_initial,&w);
778: VecDuplicate(eps->vec_initial,&xr);
779: VecDuplicate(eps->vec_initial,&xi);
780: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
782: #ifndef PETSC_USE_COMPLEX
783: if (ki == 0 ||
784: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
785: #endif
786: MatMult( A, xr, u ); /* u=A*x */
787: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
788: if (eps->isgeneralized) { MatMult( B, xr, w ); }
789: else { VecCopy( xr, w ); } /* w=B*x */
790: VecAXPY( u, -kr, w ); /* u=A*x-k*B*x */
791: }
792: VecNorm( u, NORM_2, norm);
793: #ifndef PETSC_USE_COMPLEX
794: } else {
795: MatMult( A, xr, u ); /* u=A*xr */
796: if (eps->isgeneralized) { MatMult( B, xr, v ); }
797: else { VecCopy( xr, v ); } /* v=B*xr */
798: VecAXPY( u, -kr, v ); /* u=A*xr-kr*B*xr */
799: if (eps->isgeneralized) { MatMult( B, xi, w ); }
800: else { VecCopy( xi, w ); } /* w=B*xi */
801: VecAXPY( u, ki, w ); /* u=A*xr-kr*B*xr+ki*B*xi */
802: VecNorm( u, NORM_2, &nr );
803: MatMult( A, xi, u ); /* u=A*xi */
804: VecAXPY( u, -kr, w ); /* u=A*xi-kr*B*xi */
805: VecAXPY( u, -ki, v ); /* u=A*xi-kr*B*xi-ki*B*xr */
806: VecNorm( u, NORM_2, &ni );
807: *norm = SlepcAbsEigenvalue( nr, ni );
808: }
809: #endif
811: VecDestroy(w);
812: VecDestroy(v);
813: VecDestroy(u);
814: VecDestroy(xr);
815: VecDestroy(xi);
816: return(0);
817: }
821: /*@
822: EPSComputeResidualNormLeft - Computes the norm of the residual vector associated with
823: the i-th computed left eigenvector (only available in two-sided eigensolvers).
825: Collective on EPS
827: Input Parameter:
828: . eps - the eigensolver context
829: . i - the solution index
831: Output Parameter:
832: . norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the
833: eigenvalue and y is the left eigenvector.
834: If k=0 then the residual norm is computed as ||y'A||_2.
836: Notes:
837: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
838: Eigenpairs are indexed according to the ordering criterion established
839: with EPSSetWhichEigenpairs().
841: Level: beginner
843: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
844: @*/
845: PetscErrorCode EPSComputeResidualNormLeft(EPS eps, int i, PetscReal *norm)
846: {
848: Vec u, v, w, xr, xi;
849: Mat A, B;
850: PetscScalar kr, ki;
851: #ifndef PETSC_USE_COMPLEX
852: PetscReal ni, nr;
853: #endif
854:
857: STGetOperators(eps->OP,&A,&B);
858: VecDuplicate(eps->vec_initial_left,&u);
859: VecDuplicate(eps->vec_initial_left,&v);
860: VecDuplicate(eps->vec_initial_left,&w);
861: VecDuplicate(eps->vec_initial_left,&xr);
862: VecDuplicate(eps->vec_initial_left,&xi);
863: EPSGetValue(eps,i,&kr,&ki);
864: EPSGetLeftVector(eps,i,xr,xi);
866: #ifndef PETSC_USE_COMPLEX
867: if (ki == 0 ||
868: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
869: #endif
870: MatMultTranspose( A, xr, u ); /* u=A'*x */
871: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
872: if (eps->isgeneralized) { MatMultTranspose( B, xr, w ); }
873: else { VecCopy( xr, w ); } /* w=B'*x */
874: VecAXPY( u, -kr, w); /* u=A'*x-k*B'*x */
875: }
876: VecNorm( u, NORM_2, norm);
877: #ifndef PETSC_USE_COMPLEX
878: } else {
879: MatMultTranspose( A, xr, u ); /* u=A'*xr */
880: if (eps->isgeneralized) { MatMultTranspose( B, xr, v ); }
881: else { VecCopy( xr, v ); } /* v=B'*xr */
882: VecAXPY( u, -kr, v ); /* u=A'*xr-kr*B'*xr */
883: if (eps->isgeneralized) { MatMultTranspose( B, xi, w ); }
884: else { VecCopy( xi, w ); } /* w=B'*xi */
885: VecAXPY( u, ki, w ); /* u=A'*xr-kr*B'*xr+ki*B'*xi */
886: VecNorm( u, NORM_2, &nr );
887: MatMultTranspose( A, xi, u ); /* u=A'*xi */
888: VecAXPY( u, -kr, w ); /* u=A'*xi-kr*B'*xi */
889: VecAXPY( u, -ki, v ); /* u=A'*xi-kr*B'*xi-ki*B'*xr */
890: VecNorm( u, NORM_2, &ni );
891: *norm = SlepcAbsEigenvalue( nr, ni );
892: }
893: #endif
895: VecDestroy(w);
896: VecDestroy(v);
897: VecDestroy(u);
898: VecDestroy(xr);
899: VecDestroy(xi);
900: return(0);
901: }
905: /*@
906: EPSComputeRelativeError - Computes the relative error bound associated
907: with the i-th computed eigenpair.
909: Collective on EPS
911: Input Parameter:
912: . eps - the eigensolver context
913: . i - the solution index
915: Output Parameter:
916: . error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where
917: k is the eigenvalue and x is the eigenvector.
918: If k=0 the relative error is computed as ||Ax||_2/||x||_2.
920: Level: beginner
922: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
923: @*/
924: PetscErrorCode EPSComputeRelativeError(EPS eps, int i, PetscReal *error)
925: {
927: Vec xr, xi;
928: PetscScalar kr, ki;
929: PetscReal norm, er;
930: #ifndef PETSC_USE_COMPLEX
931: Vec u;
932: PetscReal ei;
933: #endif
934:
937: EPSComputeResidualNorm(eps,i,&norm);
938: VecDuplicate(eps->vec_initial,&xr);
939: VecDuplicate(eps->vec_initial,&xi);
940: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
942: #ifndef PETSC_USE_COMPLEX
943: if (ki == 0 ||
944: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
945: #endif
946: VecNorm(xr, NORM_2, &er);
947: if (PetscAbsScalar(kr) > norm) {
948: *error = norm / (PetscAbsScalar(kr) * er);
949: } else {
950: *error = norm / er;
951: }
952: #ifndef PETSC_USE_COMPLEX
953: } else {
954: if (SlepcAbsEigenvalue(kr,ki) > norm) {
955: VecDuplicate(xi, &u);
956: VecCopy(xi, u);
957: VecAXPBY(u, kr, -ki, xr);
958: VecNorm(u, NORM_2, &er);
959: VecAXPBY(xi, kr, ki, xr);
960: VecNorm(xi, NORM_2, &ei);
961: VecDestroy(u);
962: } else {
963: VecDot(xr, xr, &er);
964: VecDot(xi, xi, &ei);
965: }
966: *error = norm / SlepcAbsEigenvalue(er, ei);
967: }
968: #endif
969:
970: VecDestroy(xr);
971: VecDestroy(xi);
972: return(0);
973: }
977: /*@
978: EPSComputeRelativeErrorLeft - Computes the relative error bound associated
979: with the i-th computed eigenvalue and left eigenvector (only available in
980: two-sided eigensolvers).
982: Collective on EPS
984: Input Parameter:
985: . eps - the eigensolver context
986: . i - the solution index
988: Output Parameter:
989: . error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where
990: k is the eigenvalue and y is the left eigenvector.
991: If k=0 the relative error is computed as ||y'A||_2/||y||_2.
993: Level: beginner
995: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
996: @*/
997: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps, int i, PetscReal *error)
998: {
1000: Vec xr, xi;
1001: PetscScalar kr, ki;
1002: PetscReal norm, er;
1003: #ifndef PETSC_USE_COMPLEX
1004: Vec u;
1005: PetscReal ei;
1006: #endif
1007:
1010: EPSComputeResidualNormLeft(eps,i,&norm);
1011: VecDuplicate(eps->vec_initial_left,&xr);
1012: VecDuplicate(eps->vec_initial_left,&xi);
1013: EPSGetValue(eps,i,&kr,&ki);
1014: EPSGetLeftVector(eps,i,xr,xi);
1016: #ifndef PETSC_USE_COMPLEX
1017: if (ki == 0 ||
1018: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1019: #endif
1020: VecNorm(xr, NORM_2, &er);
1021: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1022: *error = norm / (PetscAbsScalar(kr) * er);
1023: } else {
1024: *error = norm / er;
1025: }
1026: #ifndef PETSC_USE_COMPLEX
1027: } else {
1028: VecDuplicate(xi, &u);
1029: VecCopy(xi, u);
1030: VecAXPBY(u, kr, -ki, xr);
1031: VecNorm(u, NORM_2, &er);
1032: VecAXPBY(xi, kr, ki, xr);
1033: VecNorm(xi, NORM_2, &ei);
1034: VecDestroy(u);
1035: *error = norm / SlepcAbsEigenvalue(er, ei);
1036: }
1037: #endif
1038:
1039: VecDestroy(xr);
1040: VecDestroy(xi);
1041: return(0);
1042: }
1044: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1048: /*@
1049: EPSSortEigenvalues - Sorts a list of eigenvalues according to a certain
1050: criterion.
1052: Not Collective
1054: Input Parameters:
1055: + n - number of eigenvalue in the list
1056: . eig - pointer to the array containing the eigenvalues
1057: . eigi - imaginary part of the eigenvalues (only when using real numbers)
1058: . which - sorting criterion
1059: - nev - number of wanted eigenvalues
1061: Output Parameter:
1062: . permout - resulting permutation
1064: Notes:
1065: The result is a list of indices in the original eigenvalue array
1066: corresponding to the first nev eigenvalues sorted in the specified
1067: criterion
1069: Level: developer
1071: .seealso: EPSDenseNHEPSorted(), EPSSetWhichEigenpairs()
1072: @*/
1073: PetscErrorCode EPSSortEigenvalues(int n,PetscScalar *eig,PetscScalar *eigi,EPSWhich which,int nev,int *permout)
1074: {
1076: int i;
1077: PetscInt *perm;
1078: PetscReal *values;
1081: PetscMalloc(n*sizeof(PetscInt),&perm);
1082: PetscMalloc(n*sizeof(PetscReal),&values);
1083: for (i=0; i<n; i++) { perm[i] = i;}
1085: switch(which) {
1086: case EPS_LARGEST_MAGNITUDE:
1087: case EPS_SMALLEST_MAGNITUDE:
1088: for (i=0; i<n; i++) { values[i] = SlepcAbsEigenvalue(eig[i],eigi[i]); }
1089: break;
1090: case EPS_LARGEST_REAL:
1091: case EPS_SMALLEST_REAL:
1092: for (i=0; i<n; i++) { values[i] = PetscRealPart(eig[i]); }
1093: break;
1094: case EPS_LARGEST_IMAGINARY:
1095: case EPS_SMALLEST_IMAGINARY:
1096: #if defined(PETSC_USE_COMPLEX)
1097: for (i=0; i<n; i++) { values[i] = PetscImaginaryPart(eig[i]); }
1098: #else
1099: for (i=0; i<n; i++) { values[i] = PetscAbsReal(eigi[i]); }
1100: #endif
1101: break;
1102: default: SETERRQ(1,"Wrong value of which");
1103: }
1105: PetscSortRealWithPermutation(n,values,perm);
1107: switch(which) {
1108: case EPS_LARGEST_MAGNITUDE:
1109: case EPS_LARGEST_REAL:
1110: case EPS_LARGEST_IMAGINARY:
1111: for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
1112: break;
1113: case EPS_SMALLEST_MAGNITUDE:
1114: case EPS_SMALLEST_REAL:
1115: case EPS_SMALLEST_IMAGINARY:
1116: for (i=0; i<nev; i++) { permout[i] = perm[i]; }
1117: break;
1118: default: SETERRQ(1,"Wrong value of which");
1119: }
1121: #if !defined(PETSC_USE_COMPLEX)
1122: for (i=0; i<nev-1; i++) {
1123: if (eigi[permout[i]] != 0.0) {
1124: if (eig[permout[i]] == eig[permout[i+1]] &&
1125: eigi[permout[i]] == -eigi[permout[i+1]] &&
1126: eigi[permout[i]] < 0.0) {
1127: int tmp;
1128: SWAP(permout[i], permout[i+1], tmp);
1129: }
1130: i++;
1131: }
1132: }
1133: #endif
1135: PetscFree(values);
1136: PetscFree(perm);
1137: return(0);
1138: }
1142: /*@
1143: EPSGetStartVector - Gets a vector to be used as the starting vector
1144: in an Arnoldi or Lanczos reduction.
1146: Collective on EPS and Vec
1148: Input Parameters:
1149: + eps - the eigensolver context
1150: - i - index of the Arnoldi/Lanczos step
1152: Output Parameters:
1153: + vec - the start vector
1154: - breakdown - flag indicating that a breakdown has occurred
1156: Notes:
1157: The start vector is computed from another vector: for the first step (i=0),
1158: the initial vector is used (see EPSGetInitialVector()); otherwise a random
1159: vector is created. Then this vector is forced to be in the range of OP (only
1160: for generalized definite problems) and orthonormalized with respect to all
1161: V-vectors up to i-1.
1163: The flag breakdown is set to true if either i=0 and the vector belongs to the
1164: deflation space, or i>0 and the vector is linearly dependent with respect
1165: to the V-vectors.
1167: The caller must pass a vector already allocated with dimensions conforming
1168: to the initial vector. This vector is overwritten.
1170: Level: developer
1172: .seealso: EPSGetInitialVector()
1174: @*/
1175: PetscErrorCode EPSGetStartVector(EPS eps,int i,Vec vec,PetscTruth *breakdown)
1176: {
1178: PetscReal norm;
1179: PetscTruth lindep;
1180: Vec w;
1181:
1186: /* For the first step, use the initial vector, otherwise a random one */
1187: if (i==0) {
1188: w = eps->vec_initial;
1189: } else {
1190: VecDuplicate(eps->vec_initial,&w);
1191: SlepcVecSetRandom(w);
1192: }
1194: /* Force the vector to be in the range of OP for definite generalized problems */
1195: if (eps->ispositive) {
1196: STApply(eps->OP,w,vec);
1197: } else {
1198: VecCopy(w,vec);
1199: }
1201: /* Orthonormalize the vector with respect to previous vectors */
1202: IPOrthogonalize(eps->ip,i+eps->nds,PETSC_NULL,eps->DSV,vec,PETSC_NULL,&norm,&lindep,PETSC_NULL);
1203: if (breakdown) *breakdown = lindep;
1204: else if (lindep || norm == 0.0) {
1205: if (i==0) { SETERRQ(1,"Initial vector is zero or belongs to the deflation space"); }
1206: else { SETERRQ(1,"Unable to generate more start vectors"); }
1207: }
1208:
1209: VecScale(vec,1/norm);
1211: if (i!=0) {
1212: VecDestroy(w);
1213: }
1215: return(0);
1216: }
1219: /*@
1220: EPSGetLeftStartVector - Gets a vector to be used as the starting vector
1221: in the left recurrence of a two-sided eigensolver.
1223: Collective on EPS and Vec
1225: Input Parameters:
1226: + eps - the eigensolver context
1227: - i - index of the Arnoldi/Lanczos step
1229: Output Parameter:
1230: . vec - the start vector
1232: Notes:
1233: The start vector is computed from another vector: for the first step (i=0),
1234: the left initial vector is used (see EPSGetLeftInitialVector()); otherwise
1235: a random vector is created. Then this vector is forced to be in the range
1236: of OP' and orthonormalized with respect to all W-vectors up to i-1.
1238: The caller must pass a vector already allocated with dimensions conforming
1239: to the left initial vector. This vector is overwritten.
1241: Level: developer
1243: .seealso: EPSGetLeftInitialVector()
1245: @*/
1246: PetscErrorCode EPSGetLeftStartVector(EPS eps,int i,Vec vec)
1247: {
1249: PetscTruth breakdown;
1250: PetscReal norm;
1251: Vec w;
1252:
1257: /* For the first step, use the initial vector, otherwise a random one */
1258: if (i==0) {
1259: w = eps->vec_initial_left;
1260: }
1261: else {
1262: VecDuplicate(eps->vec_initial_left,&w);
1263: SlepcVecSetRandom(w);
1264: }
1266: /* Force the vector to be in the range of OP */
1267: STApplyTranspose(eps->OP,w,vec);
1269: /* Orthonormalize the vector with respect to previous vectors */
1270: IPOrthogonalize(eps->ip,i,PETSC_NULL,eps->W,vec,PETSC_NULL,&norm,&breakdown,PETSC_NULL);
1271: if (breakdown) {
1272: if (i==0) { SETERRQ(1,"Left initial vector is zero"); }
1273: else { SETERRQ(1,"Unable to generate more left start vectors"); }
1274: }
1275: VecScale(vec,1/norm);
1277: if (i!=0) {
1278: VecDestroy(w);
1279: }
1281: return(0);
1282: }