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-2010, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include private/epsimpl.h
26: typedef struct {
27: /* old values of eps */
28: EPSWhich
29: old_which;
30: PetscErrorCode
31: (*old_which_func)(EPS,PetscScalar,PetscScalar,PetscScalar,PetscScalar,
32: PetscInt*,void*);
33: void
34: *old_which_ctx;
35: } EPSSortForSTData;
38: PetscErrorCode EPSSortForSTFunc(EPS eps, PetscScalar ar, PetscScalar ai,
39: PetscScalar br, PetscScalar bi, PetscInt *r,
40: void *ctx)
41: {
42: EPSSortForSTData *data = (EPSSortForSTData*)ctx;
43: PetscErrorCode ierr;
47: /* Back-transform the harmonic values */
48: STBackTransform(eps->OP,1,&ar,&ai);
49: STBackTransform(eps->OP,1,&br,&bi);
51: /* Compare values using the user options for the eigenpairs selection */
52: eps->which = data->old_which;
53: eps->which_func = data->old_which_func;
54: eps->which_ctx = data->old_which_ctx;
55: EPSCompareEigenvalues(eps, ar, ai, br, bi, r);
57: /* Restore the eps values */
58: eps->which = EPS_WHICH_USER;
59: eps->which_func = EPSSortForSTFunc;
60: eps->which_ctx = data;
62: return(0);
63: }
68: /*@
69: EPSSolve - Solves the eigensystem.
71: Collective on EPS
73: Input Parameter:
74: . eps - eigensolver context obtained from EPSCreate()
76: Options Database:
77: + -eps_view - print information about the solver used
78: . -eps_view_binary - save the matrices to the default binary file
79: - -eps_plot_eigs - plot computed eigenvalues
81: Level: beginner
83: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
84: @*/
85: PetscErrorCode EPSSolve(EPS eps)
86: {
88: PetscInt i;
89: PetscReal re,im;
90: PetscScalar dot;
91: PetscTruth flg,isfold;
92: PetscViewer viewer;
93: PetscDraw draw;
94: PetscDrawSP drawsp;
95: STMatMode matmode;
96: char filename[PETSC_MAX_PATH_LEN];
97: EPSSortForSTData data;
98: Mat A,B;
99: KSP ksp;
100: Vec w,x;
101: #define NUMEXTSOLV 5
102: const EPSType solvers[NUMEXTSOLV] = { EPSARPACK, EPSBLZPACK, EPSTRLAN, EPSBLOPEX, EPSPRIMME };
107: flg = PETSC_FALSE;
108: PetscOptionsGetTruth(((PetscObject)eps)->prefix,"-eps_view_binary",&flg,PETSC_NULL);
109: if (flg) {
110: STGetOperators(eps->OP,&A,&B);
111: MatView(A,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
112: if (B) MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
113: }
115: /* reset the convergence flag from the previous solves */
116: eps->reason = EPS_CONVERGED_ITERATING;
118: /* call setup */
119: if (!eps->setupcalled){ EPSSetUp(eps); }
120: STResetOperationCounters(eps->OP);
121: IPResetOperationCounters(eps->ip);
122: eps->evecsavailable = PETSC_FALSE;
123: eps->nconv = 0;
124: eps->its = 0;
125: for (i=0;i<eps->ncv;i++) eps->eigr[i]=eps->eigi[i]=eps->errest[i]=0.0;
126: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->ncv);
128: flg = PETSC_FALSE;
129: for (i=0;i<NUMEXTSOLV && !flg;i++) {
130: PetscTypeCompare((PetscObject)eps,solvers[i],&flg);
131: }
133: PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
135: if (!flg) {
136: /* temporarily change which */
137: data.old_which = eps->which;
138: data.old_which_func = eps->which_func;
139: data.old_which_ctx = eps->which_ctx;
140: eps->which = EPS_WHICH_USER;
141: eps->which_func = EPSSortForSTFunc;
142: eps->which_ctx = &data;
143: }
145: /* call solver */
146: (*eps->ops->solve)(eps);
148: if (!flg) {
149: /* restore which */
150: eps->which = data.old_which;
151: eps->which_func = data.old_which_func;
152: eps->which_ctx = data.old_which_ctx;
153: }
155: STGetMatMode(eps->OP,&matmode);
156: if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
157: /* Purify eigenvectors before reverting operator */
158: (*eps->ops->computevectors)(eps);
159: }
160: STPostSolve(eps->OP);
161: PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
163: if (!eps->reason) {
164: SETERRQ(1,"Internal error, solver returned without setting converged reason");
165: }
167: /* Map eigenvalues back to the original problem, necessary in some
168: * spectral transformations */
169: if (eps->ops->backtransform) {
170: (*eps->ops->backtransform)(eps);
171: }
173: /* Adjust left eigenvectors in generalized problems: y = B^T y */
174: if (eps->isgeneralized && eps->leftvecs) {
175: STGetOperators(eps->OP,PETSC_NULL,&B);
176: KSPCreate(((PetscObject)eps)->comm,&ksp);
177: KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
178: KSPSetFromOptions(ksp);
179: MatGetVecs(B,PETSC_NULL,&w);
180: for (i=0;i<eps->nconv;i++) {
181: VecCopy(eps->W[i],w);
182: KSPSolveTranspose(ksp,w,eps->W[i]);
183: }
184: KSPDestroy(ksp);
185: VecDestroy(w);
186: }
188: #ifndef PETSC_USE_COMPLEX
189: /* reorder conjugate eigenvalues (positive imaginary first) */
190: for (i=0; i<eps->nconv-1; i++) {
191: if (eps->eigi[i] != 0) {
192: if (eps->eigi[i] < 0) {
193: eps->eigi[i] = -eps->eigi[i];
194: eps->eigi[i+1] = -eps->eigi[i+1];
195: if (!eps->evecsavailable) {
196: /* the next correction only works with eigenvectors */
197: (*eps->ops->computevectors)(eps);
198: }
199: VecScale(eps->V[i+1],-1.0);
200: }
201: i++;
202: }
203: }
204: #endif
206: /* quick and dirty solution for FOLD: recompute eigenvalues as Rayleigh quotients */
207: PetscTypeCompare((PetscObject)eps->OP,STFOLD,&isfold);
208: if (isfold) {
209: STGetOperators(eps->OP,&A,&B);
210: MatGetVecs(A,&w,PETSC_NULL);
211: if (!eps->evecsavailable) { (*eps->ops->computevectors)(eps); }
212: for (i=0;i<eps->nconv;i++) {
213: x = eps->V[i];
214: MatMult(A,x,w);
215: VecDot(w,x,&eps->eigr[i]);
216: if (eps->isgeneralized) {
217: MatMult(B,x,w);
218: VecDot(w,x,&dot);
219: eps->eigr[i] /= dot;
220: }
221: }
222: VecDestroy(w);
223: }
225: /* sort eigenvalues according to eps->which parameter */
226: PetscFree(eps->perm);
227: if (eps->nconv > 0) {
228: PetscMalloc(sizeof(PetscInt)*eps->nconv, &eps->perm);
229: EPSSortEigenvalues(eps, eps->nconv, eps->eigr, eps->eigi, eps->perm);
230: }
232: PetscOptionsGetString(((PetscObject)eps)->prefix,"-eps_view",filename,PETSC_MAX_PATH_LEN,&flg);
233: if (flg && !PetscPreLoadingOn) {
234: PetscViewerASCIIOpen(((PetscObject)eps)->comm,filename,&viewer);
235: EPSView(eps,viewer);
236: PetscViewerDestroy(viewer);
237: }
239: flg = PETSC_FALSE;
240: PetscOptionsGetTruth(((PetscObject)eps)->prefix,"-eps_plot_eigs",&flg,PETSC_NULL);
241: if (flg) {
242: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
243: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
244: PetscViewerDrawGetDraw(viewer,0,&draw);
245: PetscDrawSPCreate(draw,1,&drawsp);
246: for( i=0; i<eps->nconv; i++ ) {
247: #if defined(PETSC_USE_COMPLEX)
248: re = PetscRealPart(eps->eigr[i]);
249: im = PetscImaginaryPart(eps->eigi[i]);
250: #else
251: re = eps->eigr[i];
252: im = eps->eigi[i];
253: #endif
254: PetscDrawSPAddPoint(drawsp,&re,&im);
255: }
256: PetscDrawSPDraw(drawsp);
257: PetscDrawSPDestroy(drawsp);
258: PetscViewerDestroy(viewer);
259: }
261: /* Remove the initial subspaces */
262: eps->nini = 0;
263: eps->ninil = 0;
265: return(0);
266: }
270: /*@
271: EPSGetIterationNumber - Gets the current iteration number. If the
272: call to EPSSolve() is complete, then it returns the number of iterations
273: carried out by the solution method.
274:
275: Not Collective
277: Input Parameter:
278: . eps - the eigensolver context
280: Output Parameter:
281: . its - number of iterations
283: Level: intermediate
285: Note:
286: During the i-th iteration this call returns i-1. If EPSSolve() is
287: complete, then parameter "its" contains either the iteration number at
288: which convergence was successfully reached, or failure was detected.
289: Call EPSGetConvergedReason() to determine if the solver converged or
290: failed and why.
292: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
293: @*/
294: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
295: {
299: *its = eps->its;
300: return(0);
301: }
305: /*@
306: EPSGetOperationCounters - Gets the total number of operator applications,
307: inner product operations and linear iterations used by the ST object
308: during the last EPSSolve() call.
310: Not Collective
312: Input Parameter:
313: . eps - EPS context
315: Output Parameter:
316: + ops - number of operator applications
317: . dots - number of inner product operations
318: - lits - number of linear iterations
320: Notes:
321: When the eigensolver algorithm invokes STApply() then a linear system
322: must be solved (except in the case of standard eigenproblems and shift
323: transformation). The number of iterations required in this solve is
324: accumulated into a counter whose value is returned by this function.
326: These counters are reset to zero at each successive call to EPSSolve().
328: Level: intermediate
330: @*/
331: PetscErrorCode EPSGetOperationCounters(EPS eps,PetscInt* ops,PetscInt* dots,PetscInt* lits)
332: {
337: STGetOperationCounters(eps->OP,ops,lits);
338: if (dots) {
339: IPGetOperationCounters(eps->ip,dots);
340: }
341: return(0);
342: }
346: /*@
347: EPSGetConverged - Gets the number of converged eigenpairs.
349: Not Collective
351: Input Parameter:
352: . eps - the eigensolver context
353:
354: Output Parameter:
355: . nconv - number of converged eigenpairs
357: Note:
358: This function should be called after EPSSolve() has finished.
360: Level: beginner
362: .seealso: EPSSetDimensions(), EPSSolve()
363: @*/
364: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
365: {
369: *nconv = eps->nconv;
370: return(0);
371: }
376: /*@C
377: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
378: stopped.
380: Not Collective
382: Input Parameter:
383: . eps - the eigensolver context
385: Output Parameter:
386: . reason - negative value indicates diverged, positive value converged
388: Possible values for reason:
389: + EPS_CONVERGED_TOL - converged up to tolerance
390: . EPS_DIVERGED_ITS - required more than its to reach convergence
391: . EPS_DIVERGED_BREAKDOWN - generic breakdown in method
392: - EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric
394: Note:
395: Can only be called after the call to EPSSolve() is complete.
397: Level: intermediate
399: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
400: @*/
401: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
402: {
406: *reason = eps->reason;
407: return(0);
408: }
412: /*@
413: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
414: subspace.
416: Not Collective, but vectors are shared by all processors that share the EPS
418: Input Parameter:
419: . eps - the eigensolver context
420:
421: Output Parameter:
422: . v - an array of vectors
424: Notes:
425: This function should be called after EPSSolve() has finished.
427: The user should provide in v an array of nconv vectors, where nconv is
428: the value returned by EPSGetConverged().
430: The first k vectors returned in v span an invariant subspace associated
431: with the first k computed eigenvalues (note that this is not true if the
432: k-th eigenvalue is complex and matrix A is real; in this case the first
433: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
434: in X for all x in X (a similar definition applies for generalized
435: eigenproblems).
437: Level: intermediate
439: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspaceLeft()
440: @*/
441: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
442: {
444: PetscInt i;
450: if (!eps->V) {
451: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
452: }
453: if (!eps->ishermitian && eps->evecsavailable) {
454: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector,EPSComputeRelativeError or EPSComputeResidualNorm");
455: }
456: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
457: for (i=0;i<eps->nconv;i++) {
458: VecPointwiseDivide(v[i],eps->V[i],eps->D);
459: VecNormalize(v[i],PETSC_NULL);
460: }
461: }
462: else {
463: for (i=0;i<eps->nconv;i++) {
464: VecCopy(eps->V[i],v[i]);
465: }
466: }
467: return(0);
468: }
472: /*@
473: EPSGetInvariantSubspaceLeft - Gets an orthonormal basis of the computed left
474: invariant subspace (only available in two-sided eigensolvers).
476: Not Collective, but vectors are shared by all processors that share the EPS
478: Input Parameter:
479: . eps - the eigensolver context
480:
481: Output Parameter:
482: . v - an array of vectors
484: Notes:
485: This function should be called after EPSSolve() has finished.
487: The user should provide in v an array of nconv vectors, where nconv is
488: the value returned by EPSGetConverged().
490: The first k vectors returned in v span a left invariant subspace associated
491: with the first k computed eigenvalues (note that this is not true if the
492: k-th eigenvalue is complex and matrix A is real; in this case the first
493: k+1 vectors should be used). A left invariant subspace Y of A satisfies y'A
494: in Y for all y in Y (a similar definition applies for generalized
495: eigenproblems).
497: Level: intermediate
499: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
500: @*/
501: PetscErrorCode EPSGetInvariantSubspaceLeft(EPS eps, Vec *v)
502: {
504: PetscInt i;
510: if (!eps->leftvecs) {
511: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
512: }
513: if (!eps->W) {
514: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
515: }
516: if (!eps->ishermitian && eps->evecsavailable) {
517: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSGetInvariantSubspaceLeft must be called before EPSGetEigenpairLeft,EPSComputeRelativeErrorLeft or EPSComputeResidualNormLeft");
518: }
519: for (i=0;i<eps->nconv;i++) {
520: VecCopy(eps->W[i],v[i]);
521: }
522: return(0);
523: }
527: /*@
528: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
529: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
531: Not Collective, but vectors are shared by all processors that share the EPS
533: Input Parameters:
534: + eps - eigensolver context
535: - i - index of the solution
537: Output Parameters:
538: + eigr - real part of eigenvalue
539: . eigi - imaginary part of eigenvalue
540: . Vr - real part of eigenvector
541: - Vi - imaginary part of eigenvector
543: Notes:
544: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
545: configured with complex scalars the eigenvalue is stored
546: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
547: set to zero).
549: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
550: Eigenpairs are indexed according to the ordering criterion established
551: with EPSSetWhichEigenpairs().
553: The 2-norm of the eigenvector is one unless the problem is generalized
554: Hermitian. In this case the eigenvector is normalized with respect to the
555: norm defined by the B matrix.
557: Level: beginner
559: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSGetEigenvectorLeft(), EPSSolve(),
560: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
561: @*/
562: PetscErrorCode EPSGetEigenpair(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
563: {
568: if (!eps->eigr || !eps->eigi || !eps->V) {
569: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
570: }
571: if (i<0 || i>=eps->nconv) {
572: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
573: }
574: EPSGetEigenvalue(eps,i,eigr,eigi);
575: EPSGetEigenvector(eps,i,Vr,Vi);
576:
577: return(0);
578: }
582: /*@
583: EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve().
585: Not Collective
587: Input Parameters:
588: + eps - eigensolver context
589: - i - index of the solution
591: Output Parameters:
592: + eigr - real part of eigenvalue
593: - eigi - imaginary part of eigenvalue
595: Notes:
596: If the eigenvalue is real, then eigi is set to zero. If PETSc is
597: configured with complex scalars the eigenvalue is stored
598: directly in eigr (eigi is set to zero).
600: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
601: Eigenpairs are indexed according to the ordering criterion established
602: with EPSSetWhichEigenpairs().
604: Level: beginner
606: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
607: EPSGetEigenpair()
608: @*/
609: PetscErrorCode EPSGetEigenvalue(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi)
610: {
611: PetscInt k;
615: if (!eps->eigr || !eps->eigi) {
616: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
617: }
618: if (i<0 || i>=eps->nconv) {
619: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
620: }
622: if (!eps->perm) k = i;
623: else k = eps->perm[i];
624: #ifdef PETSC_USE_COMPLEX
625: if (eigr) *eigr = eps->eigr[k];
626: if (eigi) *eigi = 0;
627: #else
628: if (eigr) *eigr = eps->eigr[k];
629: if (eigi) *eigi = eps->eigi[k];
630: #endif
631:
632: return(0);
633: }
637: /*@
638: EPSGetEigenvector - Gets the i-th right eigenvector as computed by EPSSolve().
640: Not Collective, but vectors are shared by all processors that share the EPS
642: Input Parameters:
643: + eps - eigensolver context
644: - i - index of the solution
646: Output Parameters:
647: + Vr - real part of eigenvector
648: - Vi - imaginary part of eigenvector
650: Notes:
651: If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
652: configured with complex scalars the eigenvector is stored
653: directly in Vr (Vi is set to zero).
655: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
656: Eigenpairs are indexed according to the ordering criterion established
657: with EPSSetWhichEigenpairs().
659: The 2-norm of the eigenvector is one unless the problem is generalized
660: Hermitian. In this case the eigenvector is normalized with respect to the
661: norm defined by the B matrix.
663: Level: beginner
665: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
666: EPSGetEigenpair(), EPSGetEigenvectorLeft()
667: @*/
668: PetscErrorCode EPSGetEigenvector(EPS eps, PetscInt i, Vec Vr, Vec Vi)
669: {
671: PetscInt k;
675: if (!eps->V) {
676: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
677: }
678: if (i<0 || i>=eps->nconv) {
679: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
680: }
681: if (!eps->evecsavailable && (Vr || Vi) ) {
682: (*eps->ops->computevectors)(eps);
683: }
685: if (!eps->perm) k = i;
686: else k = eps->perm[i];
687: #ifdef PETSC_USE_COMPLEX
688: if (Vr) { VecCopy(eps->V[k], Vr); }
689: if (Vi) { VecSet(Vi,0.0); }
690: #else
691: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
692: if (Vr) { VecCopy(eps->V[k], Vr); }
693: if (Vi) { VecCopy(eps->V[k+1], Vi); }
694: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
695: if (Vr) { VecCopy(eps->V[k-1], Vr); }
696: if (Vi) {
697: VecCopy(eps->V[k], Vi);
698: VecScale(Vi,-1.0);
699: }
700: } else { /* real eigenvalue */
701: if (Vr) { VecCopy(eps->V[k], Vr); }
702: if (Vi) { VecSet(Vi,0.0); }
703: }
704: #endif
705:
706: return(0);
707: }
711: /*@
712: EPSGetEigenvectorLeft - Gets the i-th left eigenvector as computed by EPSSolve()
713: (only available in two-sided eigensolvers).
715: Not Collective, but vectors are shared by all processors that share the EPS
717: Input Parameters:
718: + eps - eigensolver context
719: - i - index of the solution
721: Output Parameters:
722: + Wr - real part of eigenvector
723: - Wi - imaginary part of eigenvector
725: Notes:
726: If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is
727: configured with complex scalars the eigenvector is stored
728: directly in Wr (Wi is set to zero).
730: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
731: Eigenpairs are indexed according to the ordering criterion established
732: with EPSSetWhichEigenpairs().
734: Level: beginner
736: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
737: EPSGetEigenpair(), EPSGetEigenvector()
738: @*/
739: PetscErrorCode EPSGetEigenvectorLeft(EPS eps, PetscInt i, Vec Wr, Vec Wi)
740: {
742: PetscInt k;
746: if (!eps->leftvecs) {
747: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
748: }
749: if (!eps->W) {
750: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
751: }
752: if (i<0 || i>=eps->nconv) {
753: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
754: }
755: if (!eps->evecsavailable && (Wr || Wi) ) {
756: (*eps->ops->computevectors)(eps);
757: }
759: if (!eps->perm) k = i;
760: else k = eps->perm[i];
761: #ifdef PETSC_USE_COMPLEX
762: if (Wr) { VecCopy(eps->W[k], Wr); }
763: if (Wi) { VecSet(Wi,0.0); }
764: #else
765: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
766: if (Wr) { VecCopy(eps->W[k], Wr); }
767: if (Wi) { VecCopy(eps->W[k+1], Wi); }
768: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
769: if (Wr) { VecCopy(eps->W[k-1], Wr); }
770: if (Wi) {
771: VecCopy(eps->W[k], Wi);
772: VecScale(Wi,-1.0);
773: }
774: } else { /* real eigenvalue */
775: if (Wr) { VecCopy(eps->W[k], Wr); }
776: if (Wi) { VecSet(Wi,0.0); }
777: }
778: #endif
779:
780: return(0);
781: }
785: /*@
786: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
787: computed eigenpair.
789: Not Collective
791: Input Parameter:
792: + eps - eigensolver context
793: - i - index of eigenpair
795: Output Parameter:
796: . errest - the error estimate
798: Notes:
799: This is the error estimate used internally by the eigensolver. The actual
800: error bound can be computed with EPSComputeRelativeError(). See also the users
801: manual for details.
803: Level: advanced
805: .seealso: EPSComputeRelativeError()
806: @*/
807: PetscErrorCode EPSGetErrorEstimate(EPS eps, PetscInt i, PetscReal *errest)
808: {
811: if (!eps->eigr || !eps->eigi) {
812: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
813: }
814: if (i<0 || i>=eps->nconv) {
815: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
816: }
817: if (eps->perm) i = eps->perm[i];
818: if (errest) *errest = eps->errest[i];
819: return(0);
820: }
824: /*@
825: EPSGetErrorEstimateLeft - Returns the left error estimate associated to the i-th
826: computed eigenpair (only available in two-sided eigensolvers).
828: Not Collective
830: Input Parameter:
831: + eps - eigensolver context
832: - i - index of eigenpair
834: Output Parameter:
835: . errest - the left error estimate
837: Notes:
838: This is the error estimate used internally by the eigensolver. The actual
839: error bound can be computed with EPSComputeRelativeErrorLeft(). See also the users
840: manual for details.
842: Level: advanced
844: .seealso: EPSComputeRelativeErrorLeft()
845: @*/
846: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps, PetscInt i, PetscReal *errest)
847: {
850: if (!eps->eigr || !eps->eigi) {
851: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
852: }
853: if (!eps->leftvecs) {
854: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
855: }
856: if (i<0 || i>=eps->nconv) {
857: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
858: }
859: if (eps->perm) i = eps->perm[i];
860: if (errest) *errest = eps->errest_left[i];
861: return(0);
862: }
866: /*
867: EPSComputeResidualNorm_Private - Computes the norm of the residual vector
868: associated with an eigenpair.
869: */
870: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps, PetscScalar kr, PetscScalar ki, Vec xr, Vec xi, PetscReal *norm)
871: {
873: Vec u, w;
874: Mat A, B;
875: #ifndef PETSC_USE_COMPLEX
876: Vec v;
877: PetscReal ni, nr;
878: #endif
879:
881: STGetOperators(eps->OP,&A,&B);
882: VecDuplicate(eps->V[0],&u);
883: VecDuplicate(eps->V[0],&w);
884:
885: #ifndef PETSC_USE_COMPLEX
886: if (ki == 0 ||
887: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
888: #endif
889: MatMult(A,xr,u); /* u=A*x */
890: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
891: if (eps->isgeneralized) { MatMult(B,xr,w); }
892: else { VecCopy(xr,w); } /* w=B*x */
893: VecAXPY(u,-kr,w); /* u=A*x-k*B*x */
894: }
895: VecNorm(u,NORM_2,norm);
896: #ifndef PETSC_USE_COMPLEX
897: } else {
898: VecDuplicate(eps->V[0],&v);
899: MatMult(A,xr,u); /* u=A*xr */
900: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
901: if (eps->isgeneralized) { MatMult(B,xr,v); }
902: else { VecCopy(xr,v); } /* v=B*xr */
903: VecAXPY(u,-kr,v); /* u=A*xr-kr*B*xr */
904: if (eps->isgeneralized) { MatMult(B,xi,w); }
905: else { VecCopy(xi,w); } /* w=B*xi */
906: VecAXPY(u,ki,w); /* u=A*xr-kr*B*xr+ki*B*xi */
907: }
908: VecNorm(u,NORM_2,&nr);
909: MatMult(A,xi,u); /* u=A*xi */
910: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
911: VecAXPY(u,-kr,w); /* u=A*xi-kr*B*xi */
912: VecAXPY(u,-ki,v); /* u=A*xi-kr*B*xi-ki*B*xr */
913: }
914: VecNorm(u,NORM_2,&ni);
915: *norm = SlepcAbsEigenvalue(nr,ni);
916: VecDestroy(v);
917: }
918: #endif
920: VecDestroy(w);
921: VecDestroy(u);
922: return(0);
923: }
927: /*@
928: EPSComputeResidualNorm - Computes the norm of the residual vector associated with
929: the i-th computed eigenpair.
931: Collective on EPS
933: Input Parameter:
934: . eps - the eigensolver context
935: . i - the solution index
937: Output Parameter:
938: . norm - the residual norm, computed as ||Ax-kBx||_2 where k is the
939: eigenvalue and x is the eigenvector.
940: If k=0 then the residual norm is computed as ||Ax||_2.
942: Notes:
943: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
944: Eigenpairs are indexed according to the ordering criterion established
945: with EPSSetWhichEigenpairs().
947: Level: beginner
949: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
950: @*/
951: PetscErrorCode EPSComputeResidualNorm(EPS eps, PetscInt i, PetscReal *norm)
952: {
954: Vec xr, xi;
955: PetscScalar kr, ki;
956:
960: VecDuplicate(eps->V[0],&xr);
961: VecDuplicate(eps->V[0],&xi);
962: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
963: EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,norm);
964: VecDestroy(xr);
965: VecDestroy(xi);
966: return(0);
967: }
971: /*@
972: EPSComputeResidualNormLeft - Computes the norm of the residual vector associated with
973: the i-th computed left eigenvector (only available in two-sided eigensolvers).
975: Collective on EPS
977: Input Parameter:
978: . eps - the eigensolver context
979: . i - the solution index
981: Output Parameter:
982: . norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the
983: eigenvalue and y is the left eigenvector.
984: If k=0 then the residual norm is computed as ||y'A||_2.
986: Notes:
987: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
988: Eigenpairs are indexed according to the ordering criterion established
989: with EPSSetWhichEigenpairs().
991: Level: beginner
993: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
994: @*/
995: PetscErrorCode EPSComputeResidualNormLeft(EPS eps, PetscInt i, PetscReal *norm)
996: {
998: Vec u, v, w, xr, xi;
999: Mat A, B;
1000: PetscScalar kr, ki;
1001: #ifndef PETSC_USE_COMPLEX
1002: PetscReal ni, nr;
1003: #endif
1004:
1007: if (!eps->leftvecs) {
1008: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
1009: }
1010: STGetOperators(eps->OP,&A,&B);
1011: VecDuplicate(eps->W[0],&u);
1012: VecDuplicate(eps->W[0],&v);
1013: VecDuplicate(eps->W[0],&w);
1014: VecDuplicate(eps->W[0],&xr);
1015: VecDuplicate(eps->W[0],&xi);
1016: EPSGetEigenvalue(eps,i,&kr,&ki);
1017: EPSGetEigenvectorLeft(eps,i,xr,xi);
1019: #ifndef PETSC_USE_COMPLEX
1020: if (ki == 0 ||
1021: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1022: #endif
1023: MatMultTranspose( A, xr, u ); /* u=A'*x */
1024: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1025: if (eps->isgeneralized) { MatMultTranspose( B, xr, w ); }
1026: else { VecCopy( xr, w ); } /* w=B'*x */
1027: VecAXPY( u, -kr, w); /* u=A'*x-k*B'*x */
1028: }
1029: VecNorm( u, NORM_2, norm);
1030: #ifndef PETSC_USE_COMPLEX
1031: } else {
1032: MatMultTranspose( A, xr, u ); /* u=A'*xr */
1033: if (eps->isgeneralized) { MatMultTranspose( B, xr, v ); }
1034: else { VecCopy( xr, v ); } /* v=B'*xr */
1035: VecAXPY( u, -kr, v ); /* u=A'*xr-kr*B'*xr */
1036: if (eps->isgeneralized) { MatMultTranspose( B, xi, w ); }
1037: else { VecCopy( xi, w ); } /* w=B'*xi */
1038: VecAXPY( u, ki, w ); /* u=A'*xr-kr*B'*xr+ki*B'*xi */
1039: VecNorm( u, NORM_2, &nr );
1040: MatMultTranspose( A, xi, u ); /* u=A'*xi */
1041: VecAXPY( u, -kr, w ); /* u=A'*xi-kr*B'*xi */
1042: VecAXPY( u, -ki, v ); /* u=A'*xi-kr*B'*xi-ki*B'*xr */
1043: VecNorm( u, NORM_2, &ni );
1044: *norm = SlepcAbsEigenvalue( nr, ni );
1045: }
1046: #endif
1048: VecDestroy(w);
1049: VecDestroy(v);
1050: VecDestroy(u);
1051: VecDestroy(xr);
1052: VecDestroy(xi);
1053: return(0);
1054: }
1058: /*
1059: EPSComputeRelativeError_Private - Computes the relative error bound
1060: associated with an eigenpair.
1061: */
1062: PetscErrorCode EPSComputeRelativeError_Private(EPS eps, PetscScalar kr, PetscScalar ki, Vec xr, Vec xi, PetscReal *error)
1063: {
1065: PetscReal norm, er;
1066: #ifndef PETSC_USE_COMPLEX
1067: PetscReal ei;
1068: #endif
1069:
1071: EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,&norm);
1073: #ifndef PETSC_USE_COMPLEX
1074: if (ki == 0) {
1075: #endif
1076: VecNorm(xr,NORM_2,&er);
1077: #ifndef PETSC_USE_COMPLEX
1078: } else {
1079: VecNorm(xr,NORM_2,&er);
1080: VecNorm(xi,NORM_2,&ei);
1081: er = SlepcAbsEigenvalue(er,ei);
1082: }
1083: #endif
1084: (*eps->conv_func)(eps,kr,ki,norm/er,error,eps->conv_ctx);CHKERRQ(ierr)
1085: return(0);
1086: }
1090: /*@
1091: EPSComputeRelativeError - Computes the relative error bound associated
1092: with the i-th computed eigenpair.
1094: Collective on EPS
1096: Input Parameter:
1097: . eps - the eigensolver context
1098: . i - the solution index
1100: Output Parameter:
1101: . error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where
1102: k is the eigenvalue and x is the eigenvector.
1103: If k=0 the relative error is computed as ||Ax||_2/||x||_2.
1105: Level: beginner
1107: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
1108: @*/
1109: PetscErrorCode EPSComputeRelativeError(EPS eps, PetscInt i, PetscReal *error)
1110: {
1112: Vec xr, xi;
1113: PetscScalar kr, ki;
1114:
1118: VecDuplicate(eps->V[0],&xr);
1119: VecDuplicate(eps->V[0],&xi);
1120: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
1121: EPSComputeRelativeError_Private(eps,kr,ki,xr,xi,error);
1122: VecDestroy(xr);
1123: VecDestroy(xi);
1124: return(0);
1125: }
1129: /*@
1130: EPSComputeRelativeErrorLeft - Computes the relative error bound associated
1131: with the i-th computed eigenvalue and left eigenvector (only available in
1132: two-sided eigensolvers).
1134: Collective on EPS
1136: Input Parameter:
1137: . eps - the eigensolver context
1138: . i - the solution index
1140: Output Parameter:
1141: . error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where
1142: k is the eigenvalue and y is the left eigenvector.
1143: If k=0 the relative error is computed as ||y'A||_2/||y||_2.
1145: Level: beginner
1147: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
1148: @*/
1149: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps, PetscInt i, PetscReal *error)
1150: {
1152: Vec xr, xi;
1153: PetscScalar kr, ki;
1154: PetscReal norm, er;
1155: #ifndef PETSC_USE_COMPLEX
1156: Vec u;
1157: PetscReal ei;
1158: #endif
1159:
1162: EPSComputeResidualNormLeft(eps,i,&norm);
1163: VecDuplicate(eps->W[0],&xr);
1164: VecDuplicate(eps->W[0],&xi);
1165: EPSGetEigenvalue(eps,i,&kr,&ki);
1166: EPSGetEigenvectorLeft(eps,i,xr,xi);
1168: #ifndef PETSC_USE_COMPLEX
1169: if (ki == 0 ||
1170: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1171: #endif
1172: VecNorm(xr, NORM_2, &er);
1173: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1174: *error = norm / (PetscAbsScalar(kr) * er);
1175: } else {
1176: *error = norm / er;
1177: }
1178: #ifndef PETSC_USE_COMPLEX
1179: } else {
1180: VecDuplicate(xi, &u);
1181: VecCopy(xi, u);
1182: VecAXPBY(u, kr, -ki, xr);
1183: VecNorm(u, NORM_2, &er);
1184: VecAXPBY(xi, kr, ki, xr);
1185: VecNorm(xi, NORM_2, &ei);
1186: VecDestroy(u);
1187: *error = norm / SlepcAbsEigenvalue(er, ei);
1188: }
1189: #endif
1190:
1191: VecDestroy(xr);
1192: VecDestroy(xi);
1193: return(0);
1194: }
1196: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1200: /*@
1201: EPSSortEigenvalues - Sorts a list of eigenvalues according to the criterion
1202: specified via EPSSetWhichEigenpairs().
1204: Not Collective
1206: Input Parameters:
1207: + eps - the eigensolver context
1208: . n - number of eigenvalues in the list
1209: . eigr - pointer to the array containing the eigenvalues
1210: - eigi - imaginary part of the eigenvalues (only when using real numbers)
1212: Output Parameter:
1213: . perm - resulting permutation
1215: Note:
1216: The result is a list of indices in the original eigenvalue array
1217: corresponding to the first nev eigenvalues sorted in the specified
1218: criterion.
1220: Level: developer
1222: .seealso: EPSSortEigenvaluesReal(), EPSSetWhichEigenpairs()
1223: @*/
1224: PetscErrorCode EPSSortEigenvalues(EPS eps,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
1225: {
1227: PetscScalar re,im;
1228: PetscInt i,j,result,tmp;
1231: for (i=0; i<n; i++) { perm[i] = i; }
1232: /* insertion sort */
1233: for (i=n-1; i>=0; i--) {
1234: re = eigr[perm[i]];
1235: im = eigi[perm[i]];
1236: j = i + 1;
1237: #ifndef PETSC_USE_COMPLEX
1238: if (im != 0) {
1239: /* complex eigenvalue */
1240: i--;
1241: im = eigi[perm[i]];
1242: }
1243: #endif
1244: while (j<n) {
1245: EPSCompareEigenvalues(eps,re,im,eigr[perm[j]],eigi[perm[j]],&result);
1246: if (result < 0) break;
1247: #ifndef PETSC_USE_COMPLEX
1248: /* keep together every complex conjugated eigenpair */
1249: if (im == 0) {
1250: if (eigi[perm[j]] == 0) {
1251: #endif
1252: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
1253: j++;
1254: #ifndef PETSC_USE_COMPLEX
1255: } else {
1256: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
1257: j+=2;
1258: }
1259: } else {
1260: if (eigi[perm[j]] == 0) {
1261: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
1262: j++;
1263: } else {
1264: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
1265: tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
1266: j+=2;
1267: }
1268: }
1269: #endif
1270: }
1271: }
1272: return(0);
1273: }
1277: /*@
1278: EPSSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
1279: criterion (version for real eigenvalues only).
1281: Not Collective
1283: Input Parameters:
1284: + eps - the eigensolver context
1285: . n - number of eigenvalue in the list
1286: - eig - pointer to the array containing the eigenvalues (real)
1288: Output Parameter:
1289: . perm - resulting permutation
1291: Note:
1292: The result is a list of indices in the original eigenvalue array
1293: corresponding to the first nev eigenvalues sorted in the specified
1294: criterion.
1296: Level: developer
1298: .seealso: EPSSortEigenvalues(), EPSSetWhichEigenpairs(), EPSCompareEigenvalues()
1299: @*/
1300: PetscErrorCode EPSSortEigenvaluesReal(EPS eps,PetscInt n,PetscReal *eig,PetscInt *perm)
1301: {
1303: PetscScalar re;
1304: PetscInt i,j,result,tmp;
1307: for (i=0; i<n; i++) { perm[i] = i; }
1308: /* insertion sort */
1309: for (i=1; i<n; i++) {
1310: re = eig[perm[i]];
1311: j = i-1;
1312: EPSCompareEigenvalues(eps,re,0.0,eig[perm[j]],0.0,&result);
1313: while (result<=0 && j>=0) {
1314: tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1315: if (j>=0) {
1316: EPSCompareEigenvalues(eps,re,0.0,eig[perm[j]],0.0,&result);
1317: }
1318: }
1319: }
1320: return(0);
1321: }
1325: /*@
1326: EPSCompareEigenvalues - Compares two (possibly complex) eigenvalues according
1327: to a certain criterion.
1329: Not Collective
1331: Input Parameters:
1332: + eps - the eigensolver context
1333: . ar - real part of the 1st eigenvalue
1334: . ai - imaginary part of the 1st eigenvalue
1335: . br - real part of the 2nd eigenvalue
1336: - bi - imaginary part of the 2nd eigenvalue
1338: Output Parameter:
1339: . res - result of comparison
1341: Notes:
1342: The returning parameter 'res' can be:
1343: + negative - if the 1st eigenvalue is preferred to the 2st one
1344: . zero - if both eigenvalues are equally preferred
1345: - positive - if the 2st eigenvalue is preferred to the 1st one
1347: The criterion of comparison is related to the 'which' parameter set with
1348: EPSSetWhichEigenpairs().
1350: Level: developer
1352: .seealso: EPSSortEigenvalues(), EPSSetWhichEigenpairs()
1353: @*/
1354: PetscErrorCode EPSCompareEigenvalues(EPS eps,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
1355: {
1357: PetscReal a,b;
1360: switch(eps->which) {
1361: case EPS_WHICH_USER:
1362: if (!eps->which_func) SETERRQ(1,"Undefined eigenvalue comparison function");
1363: (*eps->which_func)(eps,ar,ai,br,bi,result,eps->which_ctx);
1364: a = 0.0; b = 0.0;
1365: break;
1366: case EPS_LARGEST_MAGNITUDE:
1367: case EPS_SMALLEST_MAGNITUDE:
1368: a = SlepcAbsEigenvalue(ar,ai);
1369: b = SlepcAbsEigenvalue(br,bi);
1370: break;
1371: case EPS_LARGEST_REAL:
1372: case EPS_SMALLEST_REAL:
1373: a = PetscRealPart(ar);
1374: b = PetscRealPart(br);
1375: break;
1376: case EPS_LARGEST_IMAGINARY:
1377: case EPS_SMALLEST_IMAGINARY:
1378: #if defined(PETSC_USE_COMPLEX)
1379: a = PetscImaginaryPart(ar);
1380: b = PetscImaginaryPart(br);
1381: #else
1382: a = PetscAbsReal(ai);
1383: b = PetscAbsReal(bi);
1384: #endif
1385: break;
1386: case EPS_TARGET_MAGNITUDE:
1387: /* complex target only allowed if scalartype=complex */
1388: a = SlepcAbsEigenvalue(ar-eps->target,ai);
1389: b = SlepcAbsEigenvalue(br-eps->target,bi);
1390: break;
1391: case EPS_TARGET_REAL:
1392: a = PetscAbsReal(PetscRealPart(ar-eps->target));
1393: b = PetscAbsReal(PetscRealPart(br-eps->target));
1394: break;
1395: case EPS_TARGET_IMAGINARY:
1396: #if !defined(PETSC_USE_COMPLEX)
1397: /* complex target only allowed if scalartype=complex */
1398: a = PetscAbsReal(ai);
1399: b = PetscAbsReal(bi);
1400: #else
1401: a = PetscAbsReal(PetscImaginaryPart(ar-eps->target));
1402: b = PetscAbsReal(PetscImaginaryPart(br-eps->target));
1403: #endif
1404: break;
1405: default: SETERRQ(1,"Wrong value of which");
1406: }
1407: switch(eps->which) {
1408: case EPS_WHICH_USER:
1409: break;
1410: case EPS_LARGEST_MAGNITUDE:
1411: case EPS_LARGEST_REAL:
1412: case EPS_LARGEST_IMAGINARY:
1413: if (a<b) *result = 1;
1414: else if (a>b) *result = -1;
1415: else *result = 0;
1416: break;
1417: default:
1418: if (a>b) *result = 1;
1419: else if (a<b) *result = -1;
1420: else *result = 0;
1421: }
1422: return(0);
1423: }
1427: /*@
1428: EPSGetStartVector - Gets a suitable vector to be used as the starting vector
1429: for the recurrence that builds the right subspace.
1431: Collective on EPS and Vec
1433: Input Parameters:
1434: + eps - the eigensolver context
1435: - i - iteration number
1437: Output Parameters:
1438: + vec - the start vector
1439: - breakdown - flag indicating that a breakdown has occurred
1441: Notes:
1442: The start vector is computed from another vector: for the first step (i=0),
1443: the first initial vector is used (see EPSSetInitialSpace()); otherwise a random
1444: vector is created. Then this vector is forced to be in the range of OP (only
1445: for generalized definite problems) and orthonormalized with respect to all
1446: V-vectors up to i-1.
1448: The flag breakdown is set to true if either i=0 and the vector belongs to the
1449: deflation space, or i>0 and the vector is linearly dependent with respect
1450: to the V-vectors.
1452: The caller must pass a vector already allocated with dimensions conforming
1453: to the initial vector. This vector is overwritten.
1455: Level: developer
1457: .seealso: EPSSetInitialSpace()
1458: @*/
1459: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,Vec vec,PetscTruth *breakdown)
1460: {
1462: PetscReal norm;
1463: PetscTruth lindep;
1464: Vec w;
1465:
1470: VecDuplicate(eps->V[0],&w);
1472: /* For the first step, use the first initial vector, otherwise a random one */
1473: if (i==0 && eps->nini>0) {
1474: VecCopy(eps->V[0],w);
1475: } else {
1476: SlepcVecSetRandom(w,eps->rand);
1477: }
1479: /* Force the vector to be in the range of OP for definite generalized problems */
1480: if (eps->ispositive) {
1481: STApply(eps->OP,w,vec);
1482: } else {
1483: VecCopy(w,vec);
1484: }
1486: /* Orthonormalize the vector with respect to previous vectors */
1487: IPOrthogonalize(eps->ip,eps->nds,eps->DS,i,PETSC_NULL,eps->V,vec,PETSC_NULL,&norm,&lindep);
1488: if (breakdown) *breakdown = lindep;
1489: else if (lindep || norm == 0.0) {
1490: if (i==0) { SETERRQ(1,"Initial vector is zero or belongs to the deflation space"); }
1491: else { SETERRQ(1,"Unable to generate more start vectors"); }
1492: }
1493: VecScale(vec,1.0/norm);
1495: VecDestroy(w);
1496: return(0);
1497: }
1501: /*@
1502: EPSGetStartVectorLeft - Gets a suitable vector to be used as the starting vector
1503: in the recurrence that builds the left subspace (in methods that work with two
1504: subspaces).
1506: Collective on EPS and Vec
1508: Input Parameters:
1509: + eps - the eigensolver context
1510: - i - iteration number
1512: Output Parameter:
1513: + vec - the start vector
1514: - breakdown - flag indicating that a breakdown has occurred
1516: Notes:
1517: The start vector is computed from another vector: for the first step (i=0),
1518: the first left initial vector is used (see EPSSetInitialSpaceLeft()); otherwise
1519: a random vector is created. Then this vector is forced to be in the range
1520: of OP' and orthonormalized with respect to all W-vectors up to i-1.
1522: The flag breakdown is set to true if i>0 and the vector is linearly dependent
1523: with respect to the W-vectors.
1525: The caller must pass a vector already allocated with dimensions conforming
1526: to the left initial vector. This vector is overwritten.
1528: Level: developer
1530: .seealso: EPSSetInitialSpaceLeft()
1532: @*/
1533: PetscErrorCode EPSGetStartVectorLeft(EPS eps,PetscInt i,Vec vec,PetscTruth *breakdown)
1534: {
1536: PetscReal norm;
1537: PetscTruth lindep;
1538: Vec w;
1539:
1544: VecDuplicate(eps->W[0],&w);
1546: /* For the first step, use the first initial left vector, otherwise a random one */
1547: if (i==0 && eps->ninil>0) {
1548: VecCopy(eps->W[0],w);
1549: } else {
1550: SlepcVecSetRandom(w,eps->rand);
1551: }
1553: /* Force the vector to be in the range of OP' */
1554: STApplyTranspose(eps->OP,w,vec);
1556: /* Orthonormalize the vector with respect to previous vectors */
1557: IPOrthogonalize(eps->ip,0,PETSC_NULL,i,PETSC_NULL,eps->W,vec,PETSC_NULL,&norm,&lindep);
1558: if (breakdown) *breakdown = lindep;
1559: else if (lindep || norm == 0.0) {
1560: if (i==0) { SETERRQ(1,"Left initial vector is zero"); }
1561: else { SETERRQ(1,"Unable to generate more left start vectors"); }
1562: }
1563: VecScale(vec,1/norm);
1565: VecDestroy(w);
1566: return(0);
1567: }