Actual source code: solve.c
1: /*
2: EPS routines related to the solution process.
3: */
4: #include src/eps/epsimpl.h
8: /*@
9: EPSSolve - Solves the eigensystem.
11: Collective on EPS
13: Input Parameter:
14: . eps - eigensolver context obtained from EPSCreate()
16: Options Database:
17: + -eps_view - print information about the solver used
18: . -eps_view_binary - save the matrices to the default binary file
19: - -eps_plot_eigs - plot computed eigenvalues
21: Level: beginner
23: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
24: @*/
25: PetscErrorCode EPSSolve(EPS eps)
26: {
28: int i;
29: PetscReal re,im;
30: PetscTruth flg;
31: PetscViewer viewer;
32: PetscDraw draw;
33: PetscDrawSP drawsp;
38: PetscOptionsHasName(eps->prefix,"-eps_view_binary",&flg);
39: if (flg) {
40: Mat A,B;
41: STGetOperators(eps->OP,&A,&B);
42: MatView(A,PETSC_VIEWER_BINARY_(eps->comm));
43: if (B) MatView(B,PETSC_VIEWER_BINARY_(eps->comm));
44: }
46: /* reset the convergence flag from the previous solves */
47: eps->reason = EPS_CONVERGED_ITERATING;
49: if (!eps->setupcalled){ EPSSetUp(eps); }
50: STResetNumberLinearIterations(eps->OP);
51: eps->count_reorthog = 0;
52: eps->count_breakdown = 0;
53: eps->nv = eps->ncv;
54: eps->evecsavailable = PETSC_FALSE;
56: PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
58: switch (eps->solverclass) {
59: case EPS_ONE_SIDE:
60: (*eps->ops->solve)(eps); break;
61: case EPS_TWO_SIDE:
62: if (eps->ops->solvets) {
63: (*eps->ops->solvets)(eps); break;
64: } else {
65: SETERRQ(1,"Two-sided version unavailable for this solver");
66: }
67: default:
68: SETERRQ(1,"Wrong value of eps->solverclass");
69: }
71: STPostSolve(eps->OP);
72: PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
74: if (!eps->reason) {
75: SETERRQ(1,"Internal error, solver returned without setting converged reason");
76: }
78: /* Map eigenvalues back to the original problem, necessary in some
79: * spectral transformations */
80: (*eps->ops->backtransform)(eps);
82: /* Adjust left eigenvectors in generalized problems: y = B^T y */
83: if (eps->isgeneralized && eps->solverclass == EPS_TWO_SIDE) {
84: Mat B;
85: KSP ksp;
86: Vec w;
87: STGetOperators(eps->OP,PETSC_NULL,&B);
88: KSPCreate(eps->comm,&ksp);
89: KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
90: KSPSetFromOptions(ksp);
91: MatGetVecs(B,PETSC_NULL,&w);
92: for (i=0;i<eps->nconv;i++) {
93: VecCopy(eps->W[i],w);
94: KSPSolveTranspose(ksp,w,eps->W[i]);
95: }
96: KSPDestroy(ksp);
97: VecDestroy(w);
98: }
99:
101: #ifndef PETSC_USE_COMPLEX
102: /* reorder conjugate eigenvalues (positive imaginary first) */
103: for (i=0; i<eps->nconv-1; i++) {
104: if (eps->eigi[i] != 0) {
105: if (eps->eigi[i] < 0) {
106: eps->eigi[i] = -eps->eigi[i];
107: eps->eigi[i+1] = -eps->eigi[i+1];
108: VecScale(eps->V[i+1],-1.0);
109: }
110: i++;
111: }
112: }
113: #endif
115: /* sort eigenvalues according to eps->which parameter */
116: PetscFree(eps->perm);
117: if (eps->nconv > 0) {
118: PetscMalloc(sizeof(int)*eps->nconv, &eps->perm);
119: EPSSortEigenvalues(eps->nconv, eps->eigr, eps->eigi, eps->which, eps->nconv, eps->perm);
120: }
122: PetscOptionsHasName(eps->prefix,"-eps_view",&flg);
123: if (flg && !PetscPreLoadingOn) { EPSView(eps,PETSC_VIEWER_STDOUT_WORLD); }
125: PetscOptionsHasName(eps->prefix,"-eps_plot_eigs",&flg);
126: if (flg) {
127: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
128: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
129: PetscViewerDrawGetDraw(viewer,0,&draw);
130: PetscDrawSPCreate(draw,1,&drawsp);
131: for( i=0; i<eps->nconv; i++ ) {
132: #if defined(PETSC_USE_COMPLEX)
133: re = PetscRealPart(eps->eigr[i]);
134: im = PetscImaginaryPart(eps->eigi[i]);
135: #else
136: re = eps->eigr[i];
137: im = eps->eigi[i];
138: #endif
139: PetscDrawSPAddPoint(drawsp,&re,&im);
140: }
141: PetscDrawSPDraw(drawsp);
142: PetscDrawSPDestroy(drawsp);
143: PetscViewerDestroy(viewer);
144: }
146: return(0);
147: }
151: /*@
152: EPSGetIterationNumber - Gets the current iteration number. If the
153: call to EPSSolve() is complete, then it returns the number of iterations
154: carried out by the solution method.
155:
156: Not Collective
158: Input Parameter:
159: . eps - the eigensolver context
161: Output Parameter:
162: . its - number of iterations
164: Level: intermediate
166: Notes:
167: During the i-th iteration this call returns i-1. If EPSSolve() is
168: complete, then parameter "its" contains either the iteration number at
169: which convergence was successfully reached, or failure was detected.
170: Call EPSGetConvergedReason() to determine if the solver converged or
171: failed and why.
173: @*/
174: PetscErrorCode EPSGetIterationNumber(EPS eps,int *its)
175: {
179: *its = eps->its;
180: return(0);
181: }
185: /*@
186: EPSGetNumberLinearIterations - Gets the total number of iterations
187: required by the linear solves associated to the ST object during the
188: last EPSSolve() call.
190: Not Collective
192: Input Parameter:
193: . eps - EPS context
195: Output Parameter:
196: . lits - number of linear iterations
198: Notes:
199: When the eigensolver algorithm invokes STApply() then a linear system
200: must be solved (except in the case of standard eigenproblems and shift
201: transformation). The number of iterations required in this solve is
202: accumulated into a counter whose value is returned by this function.
204: The iteration counter is reset to zero at each successive call to EPSSolve().
206: Level: intermediate
208: @*/
209: PetscErrorCode EPSGetNumberLinearIterations(EPS eps,int* lits)
210: {
214: STGetNumberLinearIterations(eps->OP, lits);
215: return(0);
216: }
220: /*@
221: EPSGetConverged - Gets the number of converged eigenpairs.
223: Not Collective
225: Input Parameter:
226: . eps - the eigensolver context
227:
228: Output Parameter:
229: . nconv - number of converged eigenpairs
231: Note:
232: This function should be called after EPSSolve() has finished.
234: Level: beginner
236: .seealso: EPSSetDimensions()
237: @*/
238: PetscErrorCode EPSGetConverged(EPS eps,int *nconv)
239: {
242: if (nconv) *nconv = eps->nconv;
243: return(0);
244: }
249: /*@C
250: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
251: stopped.
253: Not Collective
255: Input Parameter:
256: . eps - the eigensolver context
258: Output Parameter:
259: . reason - negative value indicates diverged, positive value converged
260: (see EPSConvergedReason)
262: Possible values for reason:
263: + EPS_CONVERGED_TOL - converged up to tolerance
264: . EPS_DIVERGED_ITS - required more than its to reach convergence
265: . EPS_DIVERGED_BREAKDOWN - generic breakdown in method
266: - EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric
268: Level: intermediate
270: Notes: Can only be called after the call to EPSSolve() is complete.
272: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
273: @*/
274: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
275: {
278: *reason = eps->reason;
279: return(0);
280: }
284: /*@
285: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
286: subspace.
288: Not Collective
290: Input Parameter:
291: . eps - the eigensolver context
292:
293: Output Parameter:
294: . v - an array of vectors
296: Notes:
297: This function should be called after EPSSolve() has finished.
299: The user should provide in v an array of nconv vectors, where nconv is
300: the value returned by EPSGetConverged().
302: The first k vectors returned in v span an invariant subspace associated
303: with the first k computed eigenvalues (note that this is not true if the
304: k-th eigenvalue is complex and matrix A is real; in this case the first
305: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
306: in X for all x in X (a similar definition applies for generalized
307: eigenproblems).
309: Level: intermediate
311: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetLeftInvariantSubspace()
312: @*/
313: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
314: {
316: int i;
322: if (!eps->V) {
323: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
324: }
325: for (i=0;i<eps->nconv;i++) {
326: VecCopy(eps->V[i],v[i]);
327: }
328: return(0);
329: }
333: /*@
334: EPSGetLeftInvariantSubspace - Gets an orthonormal basis of the computed left
335: invariant subspace (only available in two-sided eigensolvers).
337: Not Collective
339: Input Parameter:
340: . eps - the eigensolver context
341:
342: Output Parameter:
343: . v - an array of vectors
345: Notes:
346: This function should be called after EPSSolve() has finished.
348: The user should provide in v an array of nconv vectors, where nconv is
349: the value returned by EPSGetConverged().
351: The first k vectors returned in v span a left invariant subspace associated
352: with the first k computed eigenvalues (note that this is not true if the
353: k-th eigenvalue is complex and matrix A is real; in this case the first
354: k+1 vectors should be used). A left invariant subspace Y of A satisfies y'A
355: in Y for all y in Y (a similar definition applies for generalized
356: eigenproblems).
358: Level: intermediate
360: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
361: @*/
362: PetscErrorCode EPSGetLeftInvariantSubspace(EPS eps, Vec *v)
363: {
365: int i;
371: if (!eps->W) {
372: if (eps->solverclass!=EPS_TWO_SIDE) {
373: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
374: } else {
375: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
376: }
377: }
378: for (i=0;i<eps->nconv;i++) {
379: VecCopy(eps->W[i],v[i]);
380: }
381: return(0);
382: }
386: /*@
387: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
388: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
390: Not Collective
392: Input Parameters:
393: + eps - eigensolver context
394: - i - index of the solution
396: Output Parameters:
397: + eigr - real part of eigenvalue
398: . eigi - imaginary part of eigenvalue
399: . Vr - real part of eigenvector
400: - Vi - imaginary part of eigenvector
402: Notes:
403: If the eigenvalue is real, then eigi and Vi are set to zero. In the
404: complex case (e.g. with BOPT=O_complex) the eigenvalue is stored
405: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
406: set to zero).
408: The index i should be a value between 0 and nconv (see EPSGetConverged()).
409: Eigenpairs are indexed according to the ordering criterion established
410: with EPSSetWhichEigenpairs().
412: Level: beginner
414: .seealso: EPSGetValue(), EPSGetRightVector(), EPSGetLeftVector(), EPSSolve(),
415: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
416: @*/
417: PetscErrorCode EPSGetEigenpair(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
418: {
423: if (!eps->eigr || !eps->eigi || !eps->V) {
424: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
425: }
426: if (i<0 || i>=eps->nconv) {
427: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
428: }
429: EPSGetValue(eps,i,eigr,eigi);
430: EPSGetRightVector(eps,i,Vr,Vi);
431:
432: return(0);
433: }
437: /*@
438: EPSGetValue - Gets the i-th eigenvalue as computed by EPSSolve().
440: Not Collective
442: Input Parameters:
443: + eps - eigensolver context
444: - i - index of the solution
446: Output Parameters:
447: + eigr - real part of eigenvalue
448: - eigi - imaginary part of eigenvalue
450: Notes:
451: If the eigenvalue is real, then eigi is set to zero. In the
452: complex case (e.g. with BOPT=O_complex) the eigenvalue is stored
453: directly in eigr (eigi is set to zero).
455: The index i should be a value between 0 and nconv (see EPSGetConverged()).
456: Eigenpairs are indexed according to the ordering criterion established
457: with EPSSetWhichEigenpairs().
459: Level: beginner
461: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
462: EPSGetEigenpair()
463: @*/
464: PetscErrorCode EPSGetValue(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi)
465: {
466: int k;
470: if (!eps->eigr || !eps->eigi) {
471: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
472: }
473: if (i<0 || i>=eps->nconv) {
474: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
475: }
477: if (!eps->perm) k = i;
478: else k = eps->perm[i];
479: #ifdef PETSC_USE_COMPLEX
480: if (eigr) *eigr = eps->eigr[k];
481: if (eigi) *eigi = 0;
482: #else
483: if (eigr) *eigr = eps->eigr[k];
484: if (eigi) *eigi = eps->eigi[k];
485: #endif
486:
487: return(0);
488: }
492: /*@
493: EPSGetRightVector - Gets the i-th right eigenvector as computed by EPSSolve().
495: Not Collective
497: Input Parameters:
498: + eps - eigensolver context
499: - i - index of the solution
501: Output Parameters:
502: + Vr - real part of eigenvector
503: - Vi - imaginary part of eigenvector
505: Notes:
506: If the corresponding eigenvalue is real, then Vi is set to zero. In the
507: complex case (e.g. with BOPT=O_complex) the eigenvector is stored
508: directly in Vr (Vi is set to zero).
510: The index i should be a value between 0 and nconv (see EPSGetConverged()).
511: Eigenpairs are indexed according to the ordering criterion established
512: with EPSSetWhichEigenpairs().
514: Level: beginner
516: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
517: EPSGetEigenpair(), EPSGetLeftVector()
518: @*/
519: PetscErrorCode EPSGetRightVector(EPS eps, int i, Vec Vr, Vec Vi)
520: {
522: int k;
526: if (!eps->V) {
527: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
528: }
529: if (i<0 || i>=eps->nconv) {
530: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
531: }
532: if (!eps->evecsavailable && (Vr || Vi) ) {
533: (*eps->ops->computevectors)(eps);
534: }
536: if (!eps->perm) k = i;
537: else k = eps->perm[i];
538: #ifdef PETSC_USE_COMPLEX
539: if (Vr) { VecCopy(eps->AV[k], Vr); }
540: if (Vi) { VecSet(Vi,0.0); }
541: #else
542: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
543: if (Vr) { VecCopy(eps->AV[k], Vr); }
544: if (Vi) { VecCopy(eps->AV[k+1], Vi); }
545: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
546: if (Vr) { VecCopy(eps->AV[k-1], Vr); }
547: if (Vi) {
548: VecCopy(eps->AV[k], Vi);
549: VecScale(Vi,-1.0);
550: }
551: } else { /* real eigenvalue */
552: if (Vr) { VecCopy(eps->AV[k], Vr); }
553: if (Vi) { VecSet(Vi,0.0); }
554: }
555: #endif
556:
557: return(0);
558: }
562: /*@
563: EPSGetLeftVector - Gets the i-th left eigenvector as computed by EPSSolve()
564: (only available in two-sided eigensolvers).
566: Not Collective
568: Input Parameters:
569: + eps - eigensolver context
570: - i - index of the solution
572: Output Parameters:
573: + Wr - real part of eigenvector
574: - Wi - imaginary part of eigenvector
576: Notes:
577: If the corresponding eigenvalue is real, then Wi is set to zero. In the
578: complex case (e.g. with BOPT=O_complex) the eigenvector is stored
579: directly in Wr (Wi is set to zero).
581: The index i should be a value between 0 and nconv (see EPSGetConverged()).
582: Eigenpairs are indexed according to the ordering criterion established
583: with EPSSetWhichEigenpairs().
585: Level: beginner
587: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
588: EPSGetEigenpair(), EPSGetLeftVector()
589: @*/
590: PetscErrorCode EPSGetLeftVector(EPS eps, int i, Vec Wr, Vec Wi)
591: {
593: int k;
597: if (!eps->W) {
598: if (eps->solverclass!=EPS_TWO_SIDE) {
599: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
600: } else {
601: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
602: }
603: }
604: if (i<0 || i>=eps->nconv) {
605: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
606: }
607: if (!eps->evecsavailable && (Wr || Wi) ) {
608: (*eps->ops->computevectors)(eps);
609: }
611: if (!eps->perm) k = i;
612: else k = eps->perm[i];
613: #ifdef PETSC_USE_COMPLEX
614: if (Wr) { VecCopy(eps->AW[k], Wr); }
615: if (Wi) { VecSet(Wi,0.0); }
616: #else
617: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
618: if (Wr) { VecCopy(eps->AW[k], Wr); }
619: if (Wi) { VecCopy(eps->AW[k+1], Wi); }
620: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
621: if (Wr) { VecCopy(eps->AW[k-1], Wr); }
622: if (Wi) {
623: VecCopy(eps->AW[k], Wi);
624: VecScale(Wi,-1.0);
625: }
626: } else { /* real eigenvalue */
627: if (Wr) { VecCopy(eps->AW[k], Wr); }
628: if (Wi) { VecSet(Wi,0.0); }
629: }
630: #endif
631:
632: return(0);
633: }
637: /*@
638: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
639: computed eigenpair.
641: Not Collective
643: Input Parameter:
644: + eps - eigensolver context
645: - i - index of eigenpair
647: Output Parameter:
648: . errest - the error estimate
650: Notes:
651: This is the error estimate used internally by the eigensolver. The actual
652: error bound can be computed with EPSComputeRelativeError(). See also the user's
653: manual for details.
655: Level: advanced
657: .seealso: EPSComputeRelativeError()
658: @*/
659: PetscErrorCode EPSGetErrorEstimate(EPS eps, int i, PetscReal *errest)
660: {
663: if (!eps->eigr || !eps->eigi) {
664: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
665: }
666: if (i<0 || i>=eps->nconv) {
667: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
668: }
669: if (eps->perm) i = eps->perm[i];
670: if (errest) *errest = eps->errest[i];
671: return(0);
672: }
676: /*@
677: EPSGetErrorEstimateLeft - Returns the left error estimate associated to the i-th
678: computed eigenpair (only available in two-sided eigensolvers).
680: Not Collective
682: Input Parameter:
683: + eps - eigensolver context
684: - i - index of eigenpair
686: Output Parameter:
687: . errest - the left error estimate
689: Notes:
690: This is the error estimate used internally by the eigensolver. The actual
691: error bound can be computed with EPSComputeRelativeErrorLeft(). See also the user's
692: manual for details.
694: Level: advanced
696: .seealso: EPSComputeRelativeErrorLeft()
697: @*/
698: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps, int i, PetscReal *errest)
699: {
702: if (!eps->eigr || !eps->eigi) {
703: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
704: }
705: if (eps->solverclass!=EPS_TWO_SIDE) {
706: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
707: }
708: if (i<0 || i>=eps->nconv) {
709: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
710: }
711: if (eps->perm) i = eps->perm[i];
712: if (errest) *errest = eps->errest_left[i];
713: return(0);
714: }
718: /*@
719: EPSComputeResidualNorm - Computes the norm of the residual vector associated with
720: the i-th computed eigenpair.
722: Collective on EPS
724: Input Parameter:
725: . eps - the eigensolver context
726: . i - the solution index
728: Output Parameter:
729: . norm - the residual norm, computed as ||Ax-kBx||_2 where k is the
730: eigenvalue and x is the eigenvector.
731: If k=0 then the residual norm is computed as ||Ax||_2.
733: Notes:
734: The index i should be a value between 0 and nconv (see EPSGetConverged()).
735: Eigenpairs are indexed according to the ordering criterion established
736: with EPSSetWhichEigenpairs().
738: Level: beginner
740: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
741: @*/
742: PetscErrorCode EPSComputeResidualNorm(EPS eps, int i, PetscReal *norm)
743: {
745: Vec u, v, w, xr, xi;
746: Mat A, B;
747: PetscScalar kr, ki;
748: #ifndef PETSC_USE_COMPLEX
749: PetscReal ni, nr;
750: #endif
751:
754: STGetOperators(eps->OP,&A,&B);
755: VecDuplicate(eps->vec_initial,&u);
756: VecDuplicate(eps->vec_initial,&v);
757: VecDuplicate(eps->vec_initial,&w);
758: VecDuplicate(eps->vec_initial,&xr);
759: VecDuplicate(eps->vec_initial,&xi);
760: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
762: #ifndef PETSC_USE_COMPLEX
763: if (ki == 0 ||
764: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
765: #endif
766: MatMult( A, xr, u ); /* u=A*x */
767: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
768: if (eps->isgeneralized) { MatMult( B, xr, w ); }
769: else { VecCopy( xr, w ); } /* w=B*x */
770: VecAXPY( u, -kr, w ); /* u=A*x-k*B*x */
771: }
772: VecNorm( u, NORM_2, norm);
773: #ifndef PETSC_USE_COMPLEX
774: } else {
775: MatMult( A, xr, u ); /* u=A*xr */
776: if (eps->isgeneralized) { MatMult( B, xr, v ); }
777: else { VecCopy( xr, v ); } /* v=B*xr */
778: VecAXPY( u, -kr, v ); /* u=A*xr-kr*B*xr */
779: if (eps->isgeneralized) { MatMult( B, xi, w ); }
780: else { VecCopy( xi, w ); } /* w=B*xi */
781: VecAXPY( u, ki, w ); /* u=A*xr-kr*B*xr+ki*B*xi */
782: VecNorm( u, NORM_2, &nr );
783: MatMult( A, xi, u ); /* u=A*xi */
784: VecAXPY( u, -kr, w ); /* u=A*xi-kr*B*xi */
785: VecAXPY( u, -ki, v ); /* u=A*xi-kr*B*xi-ki*B*xr */
786: VecNorm( u, NORM_2, &ni );
787: *norm = SlepcAbsEigenvalue( nr, ni );
788: }
789: #endif
791: VecDestroy(w);
792: VecDestroy(v);
793: VecDestroy(u);
794: VecDestroy(xr);
795: VecDestroy(xi);
796: return(0);
797: }
801: /*@
802: EPSComputeResidualNormLeft - Computes the norm of the residual vector associated with
803: the i-th computed left eigenvector (only available in two-sided eigensolvers).
805: Collective on EPS
807: Input Parameter:
808: . eps - the eigensolver context
809: . i - the solution index
811: Output Parameter:
812: . norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the
813: eigenvalue and y is the left eigenvector.
814: If k=0 then the residual norm is computed as ||y'A||_2.
816: Notes:
817: The index i should be a value between 0 and nconv (see EPSGetConverged()).
818: Eigenpairs are indexed according to the ordering criterion established
819: with EPSSetWhichEigenpairs().
821: Level: beginner
823: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
824: @*/
825: PetscErrorCode EPSComputeResidualNormLeft(EPS eps, int i, PetscReal *norm)
826: {
828: Vec u, v, w, xr, xi;
829: Mat A, B;
830: PetscScalar kr, ki;
831: #ifndef PETSC_USE_COMPLEX
832: PetscReal ni, nr;
833: #endif
834:
837: STGetOperators(eps->OP,&A,&B);
838: VecDuplicate(eps->vec_initial_left,&u);
839: VecDuplicate(eps->vec_initial_left,&v);
840: VecDuplicate(eps->vec_initial_left,&w);
841: VecDuplicate(eps->vec_initial_left,&xr);
842: VecDuplicate(eps->vec_initial_left,&xi);
843: EPSGetValue(eps,i,&kr,&ki);
844: EPSGetLeftVector(eps,i,xr,xi);
846: #ifndef PETSC_USE_COMPLEX
847: if (ki == 0 ||
848: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
849: #endif
850: MatMultTranspose( A, xr, u ); /* u=A'*x */
851: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
852: if (eps->isgeneralized) { MatMultTranspose( B, xr, w ); }
853: else { VecCopy( xr, w ); } /* w=B'*x */
854: VecAXPY( u, -kr, w); /* u=A'*x-k*B'*x */
855: }
856: VecNorm( u, NORM_2, norm);
857: #ifndef PETSC_USE_COMPLEX
858: } else {
859: MatMultTranspose( A, xr, u ); /* u=A'*xr */
860: if (eps->isgeneralized) { MatMultTranspose( B, xr, v ); }
861: else { VecCopy( xr, v ); } /* v=B'*xr */
862: VecAXPY( u, -kr, v ); /* u=A'*xr-kr*B'*xr */
863: if (eps->isgeneralized) { MatMultTranspose( B, xi, w ); }
864: else { VecCopy( xi, w ); } /* w=B'*xi */
865: VecAXPY( u, ki, w ); /* u=A'*xr-kr*B'*xr+ki*B'*xi */
866: VecNorm( u, NORM_2, &nr );
867: MatMultTranspose( A, xi, u ); /* u=A'*xi */
868: VecAXPY( u, -kr, w ); /* u=A'*xi-kr*B'*xi */
869: VecAXPY( u, -ki, v ); /* u=A'*xi-kr*B'*xi-ki*B'*xr */
870: VecNorm( u, NORM_2, &ni );
871: *norm = SlepcAbsEigenvalue( nr, ni );
872: }
873: #endif
875: VecDestroy(w);
876: VecDestroy(v);
877: VecDestroy(u);
878: VecDestroy(xr);
879: VecDestroy(xi);
880: return(0);
881: }
885: /*@
886: EPSComputeRelativeError - Computes the relative error bound associated
887: with the i-th computed eigenpair.
889: Collective on EPS
891: Input Parameter:
892: . eps - the eigensolver context
893: . i - the solution index
895: Output Parameter:
896: . error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where
897: k is the eigenvalue and x is the eigenvector.
898: If k=0 the relative error is computed as ||Ax||_2/||x||_2.
900: Level: beginner
902: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
903: @*/
904: PetscErrorCode EPSComputeRelativeError(EPS eps, int i, PetscReal *error)
905: {
907: Vec xr, xi;
908: PetscScalar kr, ki;
909: PetscReal norm, er;
910: #ifndef PETSC_USE_COMPLEX
911: Vec u;
912: PetscReal ei;
913: #endif
914:
917: EPSComputeResidualNorm(eps,i,&norm);
918: VecDuplicate(eps->vec_initial,&xr);
919: VecDuplicate(eps->vec_initial,&xi);
920: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
922: #ifndef PETSC_USE_COMPLEX
923: if (ki == 0 ||
924: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
925: #endif
926: VecNorm(xr, NORM_2, &er);
927: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
928: *error = norm / (PetscAbsScalar(kr) * er);
929: } else {
930: *error = norm / er;
931: }
932: #ifndef PETSC_USE_COMPLEX
933: } else {
934: VecDuplicate(xi, &u);
935: VecCopy(xi, u);
936: VecAXPBY(u, kr, -ki, xr);
937: VecNorm(u, NORM_2, &er);
938: VecAXPBY(xi, kr, ki, xr);
939: VecNorm(xi, NORM_2, &ei);
940: VecDestroy(u);
941: *error = norm / SlepcAbsEigenvalue(er, ei);
942: }
943: #endif
944:
945: VecDestroy(xr);
946: VecDestroy(xi);
947: return(0);
948: }
952: /*@
953: EPSComputeRelativeErrorLeft - Computes the relative error bound associated
954: with the i-th computed eigenvalue and left eigenvector (only available in
955: two-sided eigensolvers).
957: Collective on EPS
959: Input Parameter:
960: . eps - the eigensolver context
961: . i - the solution index
963: Output Parameter:
964: . error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where
965: k is the eigenvalue and y is the left eigenvector.
966: If k=0 the relative error is computed as ||y'A||_2/||y||_2.
968: Level: beginner
970: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
971: @*/
972: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps, int i, PetscReal *error)
973: {
975: Vec xr, xi;
976: PetscScalar kr, ki;
977: PetscReal norm, er;
978: #ifndef PETSC_USE_COMPLEX
979: Vec u;
980: PetscReal ei;
981: #endif
982:
985: EPSComputeResidualNormLeft(eps,i,&norm);
986: VecDuplicate(eps->vec_initial_left,&xr);
987: VecDuplicate(eps->vec_initial_left,&xi);
988: EPSGetValue(eps,i,&kr,&ki);
989: EPSGetLeftVector(eps,i,xr,xi);
991: #ifndef PETSC_USE_COMPLEX
992: if (ki == 0 ||
993: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
994: #endif
995: VecNorm(xr, NORM_2, &er);
996: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
997: *error = norm / (PetscAbsScalar(kr) * er);
998: } else {
999: *error = norm / er;
1000: }
1001: #ifndef PETSC_USE_COMPLEX
1002: } else {
1003: VecDuplicate(xi, &u);
1004: VecCopy(xi, u);
1005: VecAXPBY(u, kr, -ki, xr);
1006: VecNorm(u, NORM_2, &er);
1007: VecAXPBY(xi, kr, ki, xr);
1008: VecNorm(xi, NORM_2, &ei);
1009: VecDestroy(u);
1010: *error = norm / SlepcAbsEigenvalue(er, ei);
1011: }
1012: #endif
1013:
1014: VecDestroy(xr);
1015: VecDestroy(xi);
1016: return(0);
1017: }
1019: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1023: /*@
1024: EPSSortEigenvalues - Sorts a list of eigenvalues according to a certain
1025: criterion.
1027: Not Collective
1029: Input Parameters:
1030: + n - number of eigenvalue in the list
1031: . eig - pointer to the array containing the eigenvalues
1032: . eigi - imaginary part of the eigenvalues (only when using real numbers)
1033: . which - sorting criterion
1034: - nev - number of wanted eigenvalues
1036: Output Parameter:
1037: . permout - resulting permutation
1039: Notes:
1040: The result is a list of indices in the original eigenvalue array
1041: corresponding to the first nev eigenvalues sorted in the specified
1042: criterion
1044: Level: developer
1046: .seealso: EPSDenseNHEPSorted(), EPSSetWhichEigenpairs()
1047: @*/
1048: PetscErrorCode EPSSortEigenvalues(int n,PetscScalar *eig,PetscScalar *eigi,EPSWhich which,int nev,int *permout)
1049: {
1051: int i;
1052: PetscInt *perm;
1053: PetscReal *values;
1056: PetscMalloc(n*sizeof(PetscInt),&perm);
1057: PetscMalloc(n*sizeof(PetscReal),&values);
1058: for (i=0; i<n; i++) { perm[i] = i;}
1060: switch(which) {
1061: case EPS_LARGEST_MAGNITUDE:
1062: case EPS_SMALLEST_MAGNITUDE:
1063: for (i=0; i<n; i++) { values[i] = SlepcAbsEigenvalue(eig[i],eigi[i]); }
1064: break;
1065: case EPS_LARGEST_REAL:
1066: case EPS_SMALLEST_REAL:
1067: for (i=0; i<n; i++) { values[i] = PetscRealPart(eig[i]); }
1068: break;
1069: case EPS_LARGEST_IMAGINARY:
1070: case EPS_SMALLEST_IMAGINARY:
1071: #if defined(PETSC_USE_COMPLEX)
1072: for (i=0; i<n; i++) { values[i] = PetscImaginaryPart(eig[i]); }
1073: #else
1074: for (i=0; i<n; i++) { values[i] = PetscAbsReal(eigi[i]); }
1075: #endif
1076: break;
1077: default: SETERRQ(1,"Wrong value of which");
1078: }
1080: PetscSortRealWithPermutation(n,values,perm);
1082: switch(which) {
1083: case EPS_LARGEST_MAGNITUDE:
1084: case EPS_LARGEST_REAL:
1085: case EPS_LARGEST_IMAGINARY:
1086: for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
1087: break;
1088: case EPS_SMALLEST_MAGNITUDE:
1089: case EPS_SMALLEST_REAL:
1090: case EPS_SMALLEST_IMAGINARY:
1091: for (i=0; i<nev; i++) { permout[i] = perm[i]; }
1092: break;
1093: default: SETERRQ(1,"Wrong value of which");
1094: }
1096: #if !defined(PETSC_USE_COMPLEX)
1097: for (i=0; i<nev-1; i++) {
1098: if (eigi[permout[i]] != 0.0) {
1099: if (eig[permout[i]] == eig[permout[i+1]] &&
1100: eigi[permout[i]] == -eigi[permout[i+1]] &&
1101: eigi[permout[i]] < 0.0) {
1102: int tmp;
1103: SWAP(permout[i], permout[i+1], tmp);
1104: }
1105: i++;
1106: }
1107: }
1108: #endif
1110: PetscFree(values);
1111: PetscFree(perm);
1112: return(0);
1113: }
1117: /*@
1118: EPSGetStartVector - Gets a vector to be used as the starting vector
1119: in an Arnoldi or Lanczos reduction.
1121: Collective on EPS and Vec
1123: Input Parameters:
1124: + eps - the eigensolver context
1125: - i - index of the Arnoldi/Lanczos step
1127: Output Parameters:
1128: + vec - the start vector
1129: - breakdown - flag indicating that a breakdown has occurred
1131: Notes:
1132: The start vector is computed from another vector: for the first step (i=0),
1133: the initial vector is used (see EPSGetInitialVector()); otherwise a random
1134: vector is created. Then this vector is forced to be in the range of OP and
1135: orthonormalized with respect to all V-vectors up to i-1.
1137: The flag breakdown is set to true if either i=0 and the vector belongs to the
1138: deflation space, or i>0 and the vector is linearly dependent with respect
1139: to the V-vectors.
1141: The caller must pass a vector already allocated with dimensions conforming
1142: to the initial vector. This vector is overwritten.
1144: Level: developer
1146: .seealso: EPSGetInitialVector()
1148: @*/
1149: PetscErrorCode EPSGetStartVector(EPS eps,int i,Vec vec,PetscTruth *breakdown)
1150: {
1152: PetscReal norm;
1153: PetscTruth lindep;
1154: Vec w;
1155:
1160: /* For the first step, use the initial vector, otherwise a random one */
1161: if (i==0) {
1162: w = eps->vec_initial;
1163: } else {
1164: VecDuplicate(eps->vec_initial,&w);
1165: SlepcVecSetRandom(w);
1166: }
1168: /* Force the vector to be in the range of OP */
1169: STApply(eps->OP,w,vec);
1171: /* Orthonormalize the vector with respect to previous vectors */
1172: EPSOrthogonalize(eps,i+eps->nds,eps->DSV,vec,PETSC_NULL,&norm,&lindep);
1173: if (breakdown) *breakdown = lindep;
1174: else if (lindep) {
1175: if (i==0) { SETERRQ(1,"Initial vector is zero or belongs to the deflation space"); }
1176: else { SETERRQ(1,"Unable to generate more start vectors"); }
1177: }
1178:
1179: VecScale(vec,1/norm);
1181: if (i!=0) {
1182: VecDestroy(w);
1183: }
1185: return(0);
1186: }
1189: /*@
1190: EPSGetLeftStartVector - Gets a vector to be used as the starting vector
1191: in the left recurrence of a two-sided eigensolver.
1193: Collective on EPS and Vec
1195: Input Parameters:
1196: + eps - the eigensolver context
1197: - i - index of the Arnoldi/Lanczos step
1199: Output Parameter:
1200: . vec - the start vector
1202: Notes:
1203: The start vector is computed from another vector: for the first step (i=0),
1204: the left initial vector is used (see EPSGetLeftInitialVector()); otherwise
1205: a random vector is created. Then this vector is forced to be in the range
1206: of OP' and orthonormalized with respect to all W-vectors up to i-1.
1208: The caller must pass a vector already allocated with dimensions conforming
1209: to the left initial vector. This vector is overwritten.
1211: Level: developer
1213: .seealso: EPSGetLeftInitialVector()
1215: @*/
1216: PetscErrorCode EPSGetLeftStartVector(EPS eps,int i,Vec vec)
1217: {
1219: PetscTruth breakdown;
1220: PetscReal norm;
1221: Vec w;
1222:
1227: /* For the first step, use the initial vector, otherwise a random one */
1228: if (i==0) {
1229: w = eps->vec_initial_left;
1230: }
1231: else {
1232: VecDuplicate(eps->vec_initial_left,&w);
1233: SlepcVecSetRandom(w);
1234: }
1236: /* Force the vector to be in the range of OP */
1237: STApplyTranspose(eps->OP,w,vec);
1239: /* Orthonormalize the vector with respect to previous vectors */
1240: EPSOrthogonalize(eps,i,eps->W,vec,PETSC_NULL,&norm,&breakdown);
1241: if (breakdown) {
1242: if (i==0) { SETERRQ(1,"Left initial vector is zero"); }
1243: else { SETERRQ(1,"Unable to generate more left start vectors"); }
1244: }
1245: VecScale(vec,1/norm);
1247: if (i!=0) {
1248: VecDestroy(w);
1249: }
1251: return(0);
1252: }