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: STResetOperationCounters(eps->OP);
51: eps->nv = eps->ncv;
52: eps->evecsavailable = PETSC_FALSE;
53: eps->nconv = 0;
54: eps->its = 0;
55: for (i=0;i<eps->ncv;i++) eps->eigr[i]=eps->eigi[i]=eps->errest[i]=0.0;
56: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->nv);
58: PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
60: switch (eps->solverclass) {
61: case EPS_ONE_SIDE:
62: (*eps->ops->solve)(eps); break;
63: case EPS_TWO_SIDE:
64: if (eps->ops->solvets) {
65: (*eps->ops->solvets)(eps); break;
66: } else {
67: SETERRQ(1,"Two-sided version unavailable for this solver");
68: }
69: default:
70: SETERRQ(1,"Wrong value of eps->solverclass");
71: }
73: STPostSolve(eps->OP);
74: PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
76: if (!eps->reason) {
77: SETERRQ(1,"Internal error, solver returned without setting converged reason");
78: }
80: /* Map eigenvalues back to the original problem, necessary in some
81: * spectral transformations */
82: (*eps->ops->backtransform)(eps);
84: /* Adjust left eigenvectors in generalized problems: y = B^T y */
85: if (eps->isgeneralized && eps->solverclass == EPS_TWO_SIDE) {
86: Mat B;
87: KSP ksp;
88: Vec w;
89: STGetOperators(eps->OP,PETSC_NULL,&B);
90: KSPCreate(eps->comm,&ksp);
91: KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
92: KSPSetFromOptions(ksp);
93: MatGetVecs(B,PETSC_NULL,&w);
94: for (i=0;i<eps->nconv;i++) {
95: VecCopy(eps->W[i],w);
96: KSPSolveTranspose(ksp,w,eps->W[i]);
97: }
98: KSPDestroy(ksp);
99: VecDestroy(w);
100: }
101:
103: #ifndef PETSC_USE_COMPLEX
104: /* reorder conjugate eigenvalues (positive imaginary first) */
105: for (i=0; i<eps->nconv-1; i++) {
106: if (eps->eigi[i] != 0) {
107: if (eps->eigi[i] < 0) {
108: eps->eigi[i] = -eps->eigi[i];
109: eps->eigi[i+1] = -eps->eigi[i+1];
110: VecScale(eps->V[i+1],-1.0);
111: }
112: i++;
113: }
114: }
115: #endif
117: /* sort eigenvalues according to eps->which parameter */
118: PetscFree(eps->perm);
119: if (eps->nconv > 0) {
120: PetscMalloc(sizeof(int)*eps->nconv, &eps->perm);
121: EPSSortEigenvalues(eps->nconv, eps->eigr, eps->eigi, eps->which, eps->nconv, eps->perm);
122: }
124: PetscOptionsHasName(eps->prefix,"-eps_view",&flg);
125: if (flg && !PetscPreLoadingOn) { EPSView(eps,PETSC_VIEWER_STDOUT_WORLD); }
127: PetscOptionsHasName(eps->prefix,"-eps_plot_eigs",&flg);
128: if (flg) {
129: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
130: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
131: PetscViewerDrawGetDraw(viewer,0,&draw);
132: PetscDrawSPCreate(draw,1,&drawsp);
133: for( i=0; i<eps->nconv; i++ ) {
134: #if defined(PETSC_USE_COMPLEX)
135: re = PetscRealPart(eps->eigr[i]);
136: im = PetscImaginaryPart(eps->eigi[i]);
137: #else
138: re = eps->eigr[i];
139: im = eps->eigi[i];
140: #endif
141: PetscDrawSPAddPoint(drawsp,&re,&im);
142: }
143: PetscDrawSPDraw(drawsp);
144: PetscDrawSPDestroy(drawsp);
145: PetscViewerDestroy(viewer);
146: }
148: return(0);
149: }
153: /*@
154: EPSGetIterationNumber - Gets the current iteration number. If the
155: call to EPSSolve() is complete, then it returns the number of iterations
156: carried out by the solution method.
157:
158: Not Collective
160: Input Parameter:
161: . eps - the eigensolver context
163: Output Parameter:
164: . its - number of iterations
166: Level: intermediate
168: Notes:
169: During the i-th iteration this call returns i-1. If EPSSolve() is
170: complete, then parameter "its" contains either the iteration number at
171: which convergence was successfully reached, or failure was detected.
172: Call EPSGetConvergedReason() to determine if the solver converged or
173: failed and why.
175: @*/
176: PetscErrorCode EPSGetIterationNumber(EPS eps,int *its)
177: {
181: *its = eps->its;
182: return(0);
183: }
187: /*@
188: EPSGetOperationCounters - Gets the total number of operator applications,
189: inner product operations and linear iterations used by the ST object
190: during the last EPSSolve() call.
192: Not Collective
194: Input Parameter:
195: . eps - EPS context
197: Output Parameter:
198: + ops - number of operator applications
199: . dots - number of inner product operations
200: - lits - number of linear iterations
202: Notes:
203: When the eigensolver algorithm invokes STApply() then a linear system
204: must be solved (except in the case of standard eigenproblems and shift
205: transformation). The number of iterations required in this solve is
206: accumulated into a counter whose value is returned by this function.
208: These counters are reset to zero at each successive call to EPSSolve().
210: Level: intermediate
212: @*/
213: PetscErrorCode EPSGetOperationCounters(EPS eps,int* ops,int* dots,int* lits)
214: {
217: STGetOperationCounters(eps->OP,ops,dots,lits);
218: return(0);
219: }
223: /*@
224: EPSGetConverged - Gets the number of converged eigenpairs.
226: Not Collective
228: Input Parameter:
229: . eps - the eigensolver context
230:
231: Output Parameter:
232: . nconv - number of converged eigenpairs
234: Note:
235: This function should be called after EPSSolve() has finished.
237: Level: beginner
239: .seealso: EPSSetDimensions()
240: @*/
241: PetscErrorCode EPSGetConverged(EPS eps,int *nconv)
242: {
246: *nconv = eps->nconv;
247: return(0);
248: }
253: /*@C
254: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
255: stopped.
257: Not Collective
259: Input Parameter:
260: . eps - the eigensolver context
262: Output Parameter:
263: . reason - negative value indicates diverged, positive value converged
264: (see EPSConvergedReason)
266: Possible values for reason:
267: + EPS_CONVERGED_TOL - converged up to tolerance
268: . EPS_DIVERGED_ITS - required more than its to reach convergence
269: . EPS_DIVERGED_BREAKDOWN - generic breakdown in method
270: - EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric
272: Level: intermediate
274: Notes: Can only be called after the call to EPSSolve() is complete.
276: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
277: @*/
278: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
279: {
283: *reason = eps->reason;
284: return(0);
285: }
289: /*@
290: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
291: subspace.
293: Not Collective
295: Input Parameter:
296: . eps - the eigensolver context
297:
298: Output Parameter:
299: . v - an array of vectors
301: Notes:
302: This function should be called after EPSSolve() has finished.
304: The user should provide in v an array of nconv vectors, where nconv is
305: the value returned by EPSGetConverged().
307: The first k vectors returned in v span an invariant subspace associated
308: with the first k computed eigenvalues (note that this is not true if the
309: k-th eigenvalue is complex and matrix A is real; in this case the first
310: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
311: in X for all x in X (a similar definition applies for generalized
312: eigenproblems).
314: Level: intermediate
316: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetLeftInvariantSubspace()
317: @*/
318: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
319: {
321: int i;
327: if (!eps->V) {
328: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
329: }
330: for (i=0;i<eps->nconv;i++) {
331: VecCopy(eps->V[i],v[i]);
332: }
333: return(0);
334: }
338: /*@
339: EPSGetLeftInvariantSubspace - Gets an orthonormal basis of the computed left
340: invariant subspace (only available in two-sided eigensolvers).
342: Not Collective
344: Input Parameter:
345: . eps - the eigensolver context
346:
347: Output Parameter:
348: . v - an array of vectors
350: Notes:
351: This function should be called after EPSSolve() has finished.
353: The user should provide in v an array of nconv vectors, where nconv is
354: the value returned by EPSGetConverged().
356: The first k vectors returned in v span a left invariant subspace associated
357: with the first k computed eigenvalues (note that this is not true if the
358: k-th eigenvalue is complex and matrix A is real; in this case the first
359: k+1 vectors should be used). A left invariant subspace Y of A satisfies y'A
360: in Y for all y in Y (a similar definition applies for generalized
361: eigenproblems).
363: Level: intermediate
365: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
366: @*/
367: PetscErrorCode EPSGetLeftInvariantSubspace(EPS eps, Vec *v)
368: {
370: int i;
376: if (!eps->W) {
377: if (eps->solverclass!=EPS_TWO_SIDE) {
378: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
379: } else {
380: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
381: }
382: }
383: for (i=0;i<eps->nconv;i++) {
384: VecCopy(eps->W[i],v[i]);
385: }
386: return(0);
387: }
391: /*@
392: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
393: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
395: Not Collective
397: Input Parameters:
398: + eps - eigensolver context
399: - i - index of the solution
401: Output Parameters:
402: + eigr - real part of eigenvalue
403: . eigi - imaginary part of eigenvalue
404: . Vr - real part of eigenvector
405: - Vi - imaginary part of eigenvector
407: Notes:
408: If the eigenvalue is real, then eigi and Vi are set to zero. In the
409: complex case (e.g. with BOPT=O_complex) the eigenvalue is stored
410: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
411: set to zero).
413: The index i should be a value between 0 and nconv (see EPSGetConverged()).
414: Eigenpairs are indexed according to the ordering criterion established
415: with EPSSetWhichEigenpairs().
417: Level: beginner
419: .seealso: EPSGetValue(), EPSGetRightVector(), EPSGetLeftVector(), EPSSolve(),
420: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
421: @*/
422: PetscErrorCode EPSGetEigenpair(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
423: {
428: if (!eps->eigr || !eps->eigi || !eps->V) {
429: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
430: }
431: if (i<0 || i>=eps->nconv) {
432: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
433: }
434: EPSGetValue(eps,i,eigr,eigi);
435: EPSGetRightVector(eps,i,Vr,Vi);
436:
437: return(0);
438: }
442: /*@
443: EPSGetValue - Gets the i-th eigenvalue as computed by EPSSolve().
445: Not Collective
447: Input Parameters:
448: + eps - eigensolver context
449: - i - index of the solution
451: Output Parameters:
452: + eigr - real part of eigenvalue
453: - eigi - imaginary part of eigenvalue
455: Notes:
456: If the eigenvalue is real, then eigi is set to zero. In the
457: complex case (e.g. with BOPT=O_complex) the eigenvalue is stored
458: directly in eigr (eigi is set to zero).
460: The index i should be a value between 0 and nconv (see EPSGetConverged()).
461: Eigenpairs are indexed according to the ordering criterion established
462: with EPSSetWhichEigenpairs().
464: Level: beginner
466: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
467: EPSGetEigenpair()
468: @*/
469: PetscErrorCode EPSGetValue(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi)
470: {
471: int k;
475: if (!eps->eigr || !eps->eigi) {
476: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
477: }
478: if (i<0 || i>=eps->nconv) {
479: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
480: }
482: if (!eps->perm) k = i;
483: else k = eps->perm[i];
484: #ifdef PETSC_USE_COMPLEX
485: if (eigr) *eigr = eps->eigr[k];
486: if (eigi) *eigi = 0;
487: #else
488: if (eigr) *eigr = eps->eigr[k];
489: if (eigi) *eigi = eps->eigi[k];
490: #endif
491:
492: return(0);
493: }
497: /*@
498: EPSGetRightVector - Gets the i-th right eigenvector as computed by EPSSolve().
500: Not Collective
502: Input Parameters:
503: + eps - eigensolver context
504: - i - index of the solution
506: Output Parameters:
507: + Vr - real part of eigenvector
508: - Vi - imaginary part of eigenvector
510: Notes:
511: If the corresponding eigenvalue is real, then Vi is set to zero. In the
512: complex case (e.g. with BOPT=O_complex) the eigenvector is stored
513: directly in Vr (Vi is set to zero).
515: The index i should be a value between 0 and nconv (see EPSGetConverged()).
516: Eigenpairs are indexed according to the ordering criterion established
517: with EPSSetWhichEigenpairs().
519: Level: beginner
521: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
522: EPSGetEigenpair(), EPSGetLeftVector()
523: @*/
524: PetscErrorCode EPSGetRightVector(EPS eps, int i, Vec Vr, Vec Vi)
525: {
527: int k;
531: if (!eps->V) {
532: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
533: }
534: if (i<0 || i>=eps->nconv) {
535: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
536: }
537: if (!eps->evecsavailable && (Vr || Vi) ) {
538: (*eps->ops->computevectors)(eps);
539: }
541: if (!eps->perm) k = i;
542: else k = eps->perm[i];
543: #ifdef PETSC_USE_COMPLEX
544: if (Vr) { VecCopy(eps->AV[k], Vr); }
545: if (Vi) { VecSet(Vi,0.0); }
546: #else
547: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
548: if (Vr) { VecCopy(eps->AV[k], Vr); }
549: if (Vi) { VecCopy(eps->AV[k+1], Vi); }
550: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
551: if (Vr) { VecCopy(eps->AV[k-1], Vr); }
552: if (Vi) {
553: VecCopy(eps->AV[k], Vi);
554: VecScale(Vi,-1.0);
555: }
556: } else { /* real eigenvalue */
557: if (Vr) { VecCopy(eps->AV[k], Vr); }
558: if (Vi) { VecSet(Vi,0.0); }
559: }
560: #endif
561:
562: return(0);
563: }
567: /*@
568: EPSGetLeftVector - Gets the i-th left eigenvector as computed by EPSSolve()
569: (only available in two-sided eigensolvers).
571: Not Collective
573: Input Parameters:
574: + eps - eigensolver context
575: - i - index of the solution
577: Output Parameters:
578: + Wr - real part of eigenvector
579: - Wi - imaginary part of eigenvector
581: Notes:
582: If the corresponding eigenvalue is real, then Wi is set to zero. In the
583: complex case (e.g. with BOPT=O_complex) the eigenvector is stored
584: directly in Wr (Wi is set to zero).
586: The index i should be a value between 0 and nconv (see EPSGetConverged()).
587: Eigenpairs are indexed according to the ordering criterion established
588: with EPSSetWhichEigenpairs().
590: Level: beginner
592: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
593: EPSGetEigenpair(), EPSGetLeftVector()
594: @*/
595: PetscErrorCode EPSGetLeftVector(EPS eps, int i, Vec Wr, Vec Wi)
596: {
598: int k;
602: if (!eps->W) {
603: if (eps->solverclass!=EPS_TWO_SIDE) {
604: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
605: } else {
606: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
607: }
608: }
609: if (i<0 || i>=eps->nconv) {
610: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
611: }
612: if (!eps->evecsavailable && (Wr || Wi) ) {
613: (*eps->ops->computevectors)(eps);
614: }
616: if (!eps->perm) k = i;
617: else k = eps->perm[i];
618: #ifdef PETSC_USE_COMPLEX
619: if (Wr) { VecCopy(eps->AW[k], Wr); }
620: if (Wi) { VecSet(Wi,0.0); }
621: #else
622: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
623: if (Wr) { VecCopy(eps->AW[k], Wr); }
624: if (Wi) { VecCopy(eps->AW[k+1], Wi); }
625: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
626: if (Wr) { VecCopy(eps->AW[k-1], Wr); }
627: if (Wi) {
628: VecCopy(eps->AW[k], Wi);
629: VecScale(Wi,-1.0);
630: }
631: } else { /* real eigenvalue */
632: if (Wr) { VecCopy(eps->AW[k], Wr); }
633: if (Wi) { VecSet(Wi,0.0); }
634: }
635: #endif
636:
637: return(0);
638: }
642: /*@
643: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
644: computed eigenpair.
646: Not Collective
648: Input Parameter:
649: + eps - eigensolver context
650: - i - index of eigenpair
652: Output Parameter:
653: . errest - the error estimate
655: Notes:
656: This is the error estimate used internally by the eigensolver. The actual
657: error bound can be computed with EPSComputeRelativeError(). See also the user's
658: manual for details.
660: Level: advanced
662: .seealso: EPSComputeRelativeError()
663: @*/
664: PetscErrorCode EPSGetErrorEstimate(EPS eps, int i, PetscReal *errest)
665: {
668: if (!eps->eigr || !eps->eigi) {
669: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
670: }
671: if (i<0 || i>=eps->nconv) {
672: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
673: }
674: if (eps->perm) i = eps->perm[i];
675: if (errest) *errest = eps->errest[i];
676: return(0);
677: }
681: /*@
682: EPSGetErrorEstimateLeft - Returns the left error estimate associated to the i-th
683: computed eigenpair (only available in two-sided eigensolvers).
685: Not Collective
687: Input Parameter:
688: + eps - eigensolver context
689: - i - index of eigenpair
691: Output Parameter:
692: . errest - the left error estimate
694: Notes:
695: This is the error estimate used internally by the eigensolver. The actual
696: error bound can be computed with EPSComputeRelativeErrorLeft(). See also the user's
697: manual for details.
699: Level: advanced
701: .seealso: EPSComputeRelativeErrorLeft()
702: @*/
703: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps, int i, PetscReal *errest)
704: {
707: if (!eps->eigr || !eps->eigi) {
708: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
709: }
710: if (eps->solverclass!=EPS_TWO_SIDE) {
711: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
712: }
713: if (i<0 || i>=eps->nconv) {
714: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
715: }
716: if (eps->perm) i = eps->perm[i];
717: if (errest) *errest = eps->errest_left[i];
718: return(0);
719: }
723: /*@
724: EPSComputeResidualNorm - Computes the norm of the residual vector associated with
725: the i-th computed eigenpair.
727: Collective on EPS
729: Input Parameter:
730: . eps - the eigensolver context
731: . i - the solution index
733: Output Parameter:
734: . norm - the residual norm, computed as ||Ax-kBx||_2 where k is the
735: eigenvalue and x is the eigenvector.
736: If k=0 then the residual norm is computed as ||Ax||_2.
738: Notes:
739: The index i should be a value between 0 and nconv (see EPSGetConverged()).
740: Eigenpairs are indexed according to the ordering criterion established
741: with EPSSetWhichEigenpairs().
743: Level: beginner
745: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
746: @*/
747: PetscErrorCode EPSComputeResidualNorm(EPS eps, int i, PetscReal *norm)
748: {
750: Vec u, v, w, xr, xi;
751: Mat A, B;
752: PetscScalar kr, ki;
753: #ifndef PETSC_USE_COMPLEX
754: PetscReal ni, nr;
755: #endif
756:
759: STGetOperators(eps->OP,&A,&B);
760: VecDuplicate(eps->vec_initial,&u);
761: VecDuplicate(eps->vec_initial,&v);
762: VecDuplicate(eps->vec_initial,&w);
763: VecDuplicate(eps->vec_initial,&xr);
764: VecDuplicate(eps->vec_initial,&xi);
765: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
767: #ifndef PETSC_USE_COMPLEX
768: if (ki == 0 ||
769: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
770: #endif
771: MatMult( A, xr, u ); /* u=A*x */
772: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
773: if (eps->isgeneralized) { MatMult( B, xr, w ); }
774: else { VecCopy( xr, w ); } /* w=B*x */
775: VecAXPY( u, -kr, w ); /* u=A*x-k*B*x */
776: }
777: VecNorm( u, NORM_2, norm);
778: #ifndef PETSC_USE_COMPLEX
779: } else {
780: MatMult( A, xr, u ); /* u=A*xr */
781: if (eps->isgeneralized) { MatMult( B, xr, v ); }
782: else { VecCopy( xr, v ); } /* v=B*xr */
783: VecAXPY( u, -kr, v ); /* u=A*xr-kr*B*xr */
784: if (eps->isgeneralized) { MatMult( B, xi, w ); }
785: else { VecCopy( xi, w ); } /* w=B*xi */
786: VecAXPY( u, ki, w ); /* u=A*xr-kr*B*xr+ki*B*xi */
787: VecNorm( u, NORM_2, &nr );
788: MatMult( A, xi, u ); /* u=A*xi */
789: VecAXPY( u, -kr, w ); /* u=A*xi-kr*B*xi */
790: VecAXPY( u, -ki, v ); /* u=A*xi-kr*B*xi-ki*B*xr */
791: VecNorm( u, NORM_2, &ni );
792: *norm = SlepcAbsEigenvalue( nr, ni );
793: }
794: #endif
796: VecDestroy(w);
797: VecDestroy(v);
798: VecDestroy(u);
799: VecDestroy(xr);
800: VecDestroy(xi);
801: return(0);
802: }
806: /*@
807: EPSComputeResidualNormLeft - Computes the norm of the residual vector associated with
808: the i-th computed left eigenvector (only available in two-sided eigensolvers).
810: Collective on EPS
812: Input Parameter:
813: . eps - the eigensolver context
814: . i - the solution index
816: Output Parameter:
817: . norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the
818: eigenvalue and y is the left eigenvector.
819: If k=0 then the residual norm is computed as ||y'A||_2.
821: Notes:
822: The index i should be a value between 0 and nconv (see EPSGetConverged()).
823: Eigenpairs are indexed according to the ordering criterion established
824: with EPSSetWhichEigenpairs().
826: Level: beginner
828: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
829: @*/
830: PetscErrorCode EPSComputeResidualNormLeft(EPS eps, int i, PetscReal *norm)
831: {
833: Vec u, v, w, xr, xi;
834: Mat A, B;
835: PetscScalar kr, ki;
836: #ifndef PETSC_USE_COMPLEX
837: PetscReal ni, nr;
838: #endif
839:
842: STGetOperators(eps->OP,&A,&B);
843: VecDuplicate(eps->vec_initial_left,&u);
844: VecDuplicate(eps->vec_initial_left,&v);
845: VecDuplicate(eps->vec_initial_left,&w);
846: VecDuplicate(eps->vec_initial_left,&xr);
847: VecDuplicate(eps->vec_initial_left,&xi);
848: EPSGetValue(eps,i,&kr,&ki);
849: EPSGetLeftVector(eps,i,xr,xi);
851: #ifndef PETSC_USE_COMPLEX
852: if (ki == 0 ||
853: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
854: #endif
855: MatMultTranspose( A, xr, u ); /* u=A'*x */
856: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
857: if (eps->isgeneralized) { MatMultTranspose( B, xr, w ); }
858: else { VecCopy( xr, w ); } /* w=B'*x */
859: VecAXPY( u, -kr, w); /* u=A'*x-k*B'*x */
860: }
861: VecNorm( u, NORM_2, norm);
862: #ifndef PETSC_USE_COMPLEX
863: } else {
864: MatMultTranspose( A, xr, u ); /* u=A'*xr */
865: if (eps->isgeneralized) { MatMultTranspose( B, xr, v ); }
866: else { VecCopy( xr, v ); } /* v=B'*xr */
867: VecAXPY( u, -kr, v ); /* u=A'*xr-kr*B'*xr */
868: if (eps->isgeneralized) { MatMultTranspose( B, xi, w ); }
869: else { VecCopy( xi, w ); } /* w=B'*xi */
870: VecAXPY( u, ki, w ); /* u=A'*xr-kr*B'*xr+ki*B'*xi */
871: VecNorm( u, NORM_2, &nr );
872: MatMultTranspose( A, xi, u ); /* u=A'*xi */
873: VecAXPY( u, -kr, w ); /* u=A'*xi-kr*B'*xi */
874: VecAXPY( u, -ki, v ); /* u=A'*xi-kr*B'*xi-ki*B'*xr */
875: VecNorm( u, NORM_2, &ni );
876: *norm = SlepcAbsEigenvalue( nr, ni );
877: }
878: #endif
880: VecDestroy(w);
881: VecDestroy(v);
882: VecDestroy(u);
883: VecDestroy(xr);
884: VecDestroy(xi);
885: return(0);
886: }
890: /*@
891: EPSComputeRelativeError - Computes the relative error bound associated
892: with the i-th computed eigenpair.
894: Collective on EPS
896: Input Parameter:
897: . eps - the eigensolver context
898: . i - the solution index
900: Output Parameter:
901: . error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where
902: k is the eigenvalue and x is the eigenvector.
903: If k=0 the relative error is computed as ||Ax||_2/||x||_2.
905: Level: beginner
907: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
908: @*/
909: PetscErrorCode EPSComputeRelativeError(EPS eps, int i, PetscReal *error)
910: {
912: Vec xr, xi;
913: PetscScalar kr, ki;
914: PetscReal norm, er;
915: #ifndef PETSC_USE_COMPLEX
916: Vec u;
917: PetscReal ei;
918: #endif
919:
922: EPSComputeResidualNorm(eps,i,&norm);
923: VecDuplicate(eps->vec_initial,&xr);
924: VecDuplicate(eps->vec_initial,&xi);
925: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
927: #ifndef PETSC_USE_COMPLEX
928: if (ki == 0 ||
929: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
930: #endif
931: VecNorm(xr, NORM_2, &er);
932: if (PetscAbsScalar(kr) > norm) {
933: *error = norm / (PetscAbsScalar(kr) * er);
934: } else {
935: *error = norm / er;
936: }
937: #ifndef PETSC_USE_COMPLEX
938: } else {
939: if (SlepcAbsEigenvalue(kr,ki) > norm) {
940: VecDuplicate(xi, &u);
941: VecCopy(xi, u);
942: VecAXPBY(u, kr, -ki, xr);
943: VecNorm(u, NORM_2, &er);
944: VecAXPBY(xi, kr, ki, xr);
945: VecNorm(xi, NORM_2, &ei);
946: VecDestroy(u);
947: } else {
948: VecDot(xr, xr, &er);
949: VecDot(xi, xi, &ei);
950: }
951: *error = norm / SlepcAbsEigenvalue(er, ei);
952: }
953: #endif
954:
955: VecDestroy(xr);
956: VecDestroy(xi);
957: return(0);
958: }
962: /*@
963: EPSComputeRelativeErrorLeft - Computes the relative error bound associated
964: with the i-th computed eigenvalue and left eigenvector (only available in
965: two-sided eigensolvers).
967: Collective on EPS
969: Input Parameter:
970: . eps - the eigensolver context
971: . i - the solution index
973: Output Parameter:
974: . error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where
975: k is the eigenvalue and y is the left eigenvector.
976: If k=0 the relative error is computed as ||y'A||_2/||y||_2.
978: Level: beginner
980: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
981: @*/
982: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps, int i, PetscReal *error)
983: {
985: Vec xr, xi;
986: PetscScalar kr, ki;
987: PetscReal norm, er;
988: #ifndef PETSC_USE_COMPLEX
989: Vec u;
990: PetscReal ei;
991: #endif
992:
995: EPSComputeResidualNormLeft(eps,i,&norm);
996: VecDuplicate(eps->vec_initial_left,&xr);
997: VecDuplicate(eps->vec_initial_left,&xi);
998: EPSGetValue(eps,i,&kr,&ki);
999: EPSGetLeftVector(eps,i,xr,xi);
1001: #ifndef PETSC_USE_COMPLEX
1002: if (ki == 0 ||
1003: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1004: #endif
1005: VecNorm(xr, NORM_2, &er);
1006: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1007: *error = norm / (PetscAbsScalar(kr) * er);
1008: } else {
1009: *error = norm / er;
1010: }
1011: #ifndef PETSC_USE_COMPLEX
1012: } else {
1013: VecDuplicate(xi, &u);
1014: VecCopy(xi, u);
1015: VecAXPBY(u, kr, -ki, xr);
1016: VecNorm(u, NORM_2, &er);
1017: VecAXPBY(xi, kr, ki, xr);
1018: VecNorm(xi, NORM_2, &ei);
1019: VecDestroy(u);
1020: *error = norm / SlepcAbsEigenvalue(er, ei);
1021: }
1022: #endif
1023:
1024: VecDestroy(xr);
1025: VecDestroy(xi);
1026: return(0);
1027: }
1029: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1033: /*@
1034: EPSSortEigenvalues - Sorts a list of eigenvalues according to a certain
1035: criterion.
1037: Not Collective
1039: Input Parameters:
1040: + n - number of eigenvalue in the list
1041: . eig - pointer to the array containing the eigenvalues
1042: . eigi - imaginary part of the eigenvalues (only when using real numbers)
1043: . which - sorting criterion
1044: - nev - number of wanted eigenvalues
1046: Output Parameter:
1047: . permout - resulting permutation
1049: Notes:
1050: The result is a list of indices in the original eigenvalue array
1051: corresponding to the first nev eigenvalues sorted in the specified
1052: criterion
1054: Level: developer
1056: .seealso: EPSDenseNHEPSorted(), EPSSetWhichEigenpairs()
1057: @*/
1058: PetscErrorCode EPSSortEigenvalues(int n,PetscScalar *eig,PetscScalar *eigi,EPSWhich which,int nev,int *permout)
1059: {
1061: int i;
1062: PetscInt *perm;
1063: PetscReal *values;
1066: PetscMalloc(n*sizeof(PetscInt),&perm);
1067: PetscMalloc(n*sizeof(PetscReal),&values);
1068: for (i=0; i<n; i++) { perm[i] = i;}
1070: switch(which) {
1071: case EPS_LARGEST_MAGNITUDE:
1072: case EPS_SMALLEST_MAGNITUDE:
1073: for (i=0; i<n; i++) { values[i] = SlepcAbsEigenvalue(eig[i],eigi[i]); }
1074: break;
1075: case EPS_LARGEST_REAL:
1076: case EPS_SMALLEST_REAL:
1077: for (i=0; i<n; i++) { values[i] = PetscRealPart(eig[i]); }
1078: break;
1079: case EPS_LARGEST_IMAGINARY:
1080: case EPS_SMALLEST_IMAGINARY:
1081: #if defined(PETSC_USE_COMPLEX)
1082: for (i=0; i<n; i++) { values[i] = PetscImaginaryPart(eig[i]); }
1083: #else
1084: for (i=0; i<n; i++) { values[i] = PetscAbsReal(eigi[i]); }
1085: #endif
1086: break;
1087: default: SETERRQ(1,"Wrong value of which");
1088: }
1090: PetscSortRealWithPermutation(n,values,perm);
1092: switch(which) {
1093: case EPS_LARGEST_MAGNITUDE:
1094: case EPS_LARGEST_REAL:
1095: case EPS_LARGEST_IMAGINARY:
1096: for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
1097: break;
1098: case EPS_SMALLEST_MAGNITUDE:
1099: case EPS_SMALLEST_REAL:
1100: case EPS_SMALLEST_IMAGINARY:
1101: for (i=0; i<nev; i++) { permout[i] = perm[i]; }
1102: break;
1103: default: SETERRQ(1,"Wrong value of which");
1104: }
1106: #if !defined(PETSC_USE_COMPLEX)
1107: for (i=0; i<nev-1; i++) {
1108: if (eigi[permout[i]] != 0.0) {
1109: if (eig[permout[i]] == eig[permout[i+1]] &&
1110: eigi[permout[i]] == -eigi[permout[i+1]] &&
1111: eigi[permout[i]] < 0.0) {
1112: int tmp;
1113: SWAP(permout[i], permout[i+1], tmp);
1114: }
1115: i++;
1116: }
1117: }
1118: #endif
1120: PetscFree(values);
1121: PetscFree(perm);
1122: return(0);
1123: }
1127: /*@
1128: EPSGetStartVector - Gets a vector to be used as the starting vector
1129: in an Arnoldi or Lanczos reduction.
1131: Collective on EPS and Vec
1133: Input Parameters:
1134: + eps - the eigensolver context
1135: - i - index of the Arnoldi/Lanczos step
1137: Output Parameters:
1138: + vec - the start vector
1139: - breakdown - flag indicating that a breakdown has occurred
1141: Notes:
1142: The start vector is computed from another vector: for the first step (i=0),
1143: the initial vector is used (see EPSGetInitialVector()); otherwise a random
1144: vector is created. Then this vector is forced to be in the range of OP and
1145: orthonormalized with respect to all V-vectors up to i-1.
1147: The flag breakdown is set to true if either i=0 and the vector belongs to the
1148: deflation space, or i>0 and the vector is linearly dependent with respect
1149: to the V-vectors.
1151: The caller must pass a vector already allocated with dimensions conforming
1152: to the initial vector. This vector is overwritten.
1154: Level: developer
1156: .seealso: EPSGetInitialVector()
1158: @*/
1159: PetscErrorCode EPSGetStartVector(EPS eps,int i,Vec vec,PetscTruth *breakdown)
1160: {
1162: PetscReal norm;
1163: PetscTruth lindep;
1164: Vec w;
1165:
1170: /* For the first step, use the initial vector, otherwise a random one */
1171: if (i==0) {
1172: w = eps->vec_initial;
1173: } else {
1174: VecDuplicate(eps->vec_initial,&w);
1175: SlepcVecSetRandom(w);
1176: }
1178: /* Force the vector to be in the range of OP */
1179: STApply(eps->OP,w,vec);
1181: /* Orthonormalize the vector with respect to previous vectors */
1182: EPSOrthogonalize(eps,i+eps->nds,PETSC_NULL,eps->DSV,vec,PETSC_NULL,&norm,&lindep);
1183: if (breakdown) *breakdown = lindep;
1184: else if (lindep || norm == 0.0) {
1185: if (i==0) { SETERRQ(1,"Initial vector is zero or belongs to the deflation space"); }
1186: else { SETERRQ(1,"Unable to generate more start vectors"); }
1187: }
1188:
1189: VecScale(vec,1/norm);
1191: if (i!=0) {
1192: VecDestroy(w);
1193: }
1195: return(0);
1196: }
1199: /*@
1200: EPSGetLeftStartVector - Gets a vector to be used as the starting vector
1201: in the left recurrence of a two-sided eigensolver.
1203: Collective on EPS and Vec
1205: Input Parameters:
1206: + eps - the eigensolver context
1207: - i - index of the Arnoldi/Lanczos step
1209: Output Parameter:
1210: . vec - the start vector
1212: Notes:
1213: The start vector is computed from another vector: for the first step (i=0),
1214: the left initial vector is used (see EPSGetLeftInitialVector()); otherwise
1215: a random vector is created. Then this vector is forced to be in the range
1216: of OP' and orthonormalized with respect to all W-vectors up to i-1.
1218: The caller must pass a vector already allocated with dimensions conforming
1219: to the left initial vector. This vector is overwritten.
1221: Level: developer
1223: .seealso: EPSGetLeftInitialVector()
1225: @*/
1226: PetscErrorCode EPSGetLeftStartVector(EPS eps,int i,Vec vec)
1227: {
1229: PetscTruth breakdown;
1230: PetscReal norm;
1231: Vec w;
1232:
1237: /* For the first step, use the initial vector, otherwise a random one */
1238: if (i==0) {
1239: w = eps->vec_initial_left;
1240: }
1241: else {
1242: VecDuplicate(eps->vec_initial_left,&w);
1243: SlepcVecSetRandom(w);
1244: }
1246: /* Force the vector to be in the range of OP */
1247: STApplyTranspose(eps->OP,w,vec);
1249: /* Orthonormalize the vector with respect to previous vectors */
1250: EPSOrthogonalize(eps,i,PETSC_NULL,eps->W,vec,PETSC_NULL,&norm,&breakdown);
1251: if (breakdown) {
1252: if (i==0) { SETERRQ(1,"Left initial vector is zero"); }
1253: else { SETERRQ(1,"Unable to generate more left start vectors"); }
1254: }
1255: VecScale(vec,1/norm);
1257: if (i!=0) {
1258: VecDestroy(w);
1259: }
1261: return(0);
1262: }