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-2012, Universitat 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 <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
26: typedef struct {
27: PetscErrorCode (*which_func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar, PetscInt*,void*);
28: void *which_ctx;
29: ST st;
30: } EPSSortForSTData;
34: PetscErrorCode EPSSortForSTFunc(PetscScalar ar,PetscScalar ai,
35: PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
36: {
37: EPSSortForSTData *data = (EPSSortForSTData*)ctx;
38: PetscErrorCode ierr;
41: STBackTransform(data->st,1,&ar,&ai);
42: STBackTransform(data->st,1,&br,&bi);
43: (*data->which_func)(ar,ai,br,bi,r,data->which_ctx);
44: return(0);
45: }
49: /*@
50: EPSSolve - Solves the eigensystem.
52: Collective on EPS
54: Input Parameter:
55: . eps - eigensolver context obtained from EPSCreate()
57: Options Database:
58: + -eps_view - print information about the solver used
59: . -eps_view_binary - save the matrices to the default binary file
60: - -eps_plot_eigs - plot computed eigenvalues
62: Level: beginner
64: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
65: @*/
66: PetscErrorCode EPSSolve(EPS eps)
67: {
69: PetscInt i;
70: PetscReal re,im;
71: PetscScalar dot;
72: PetscBool flg,isfold,iscayley,viewed=PETSC_FALSE;
73: PetscViewer viewer;
74: PetscDraw draw;
75: PetscDrawSP drawsp;
76: STMatMode matmode;
77: char filename[PETSC_MAX_PATH_LEN];
78: char view[10];
79: EPSSortForSTData data;
80: Mat A,B;
81: KSP ksp;
82: Vec w,x;
87: flg = PETSC_FALSE;
88: PetscOptionsGetBool(((PetscObject)eps)->prefix,"-eps_view_binary",&flg,PETSC_NULL);
89: if (flg) {
90: STGetOperators(eps->OP,&A,&B);
91: MatView(A,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
92: if (B) MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
93: }
95: PetscOptionsGetString(((PetscObject)eps)->prefix,"-eps_view",view,10,&flg);
96: if (flg) {
97: PetscStrcmp(view,"before",&viewed);
98: if (viewed){
99: PetscViewer viewer;
100: PetscViewerASCIIGetStdout(((PetscObject)eps)->comm,&viewer);
101: EPSView(eps,viewer);
102: }
103: }
105: /* reset the convergence flag from the previous solves */
106: eps->reason = EPS_CONVERGED_ITERATING;
108: /* call setup */
109: if (!eps->setupcalled) { EPSSetUp(eps); }
110: eps->evecsavailable = PETSC_FALSE;
111: eps->nconv = 0;
112: eps->its = 0;
113: for (i=0;i<eps->ncv;i++) eps->eigr[i]=eps->eigi[i]=eps->errest[i]=0.0;
114: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->ncv);
116: PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
118: PetscObjectTypeCompareAny((PetscObject)eps,&flg,EPSARPACK,EPSBLZPACK,EPSTRLAN,EPSBLOPEX,EPSPRIMME,"");
119: if (!flg) {
120: /* temporarily change eigenvalue comparison function */
121: data.which_func = eps->which_func;
122: data.which_ctx = eps->which_ctx;
123: data.st = eps->OP;
124: eps->which_func = (eps->which==EPS_ALL)? SlepcCompareLargestMagnitude: EPSSortForSTFunc;
125: eps->which_ctx = (eps->which==EPS_ALL)? PETSC_NULL: &data;
126: }
127: DSSetEigenvalueComparison(eps->ds,eps->which_func,eps->which_ctx);
129: /* call solver */
130: (*eps->ops->solve)(eps);
132: if (!flg) {
133: /* restore comparison function */
134: eps->which_func = data.which_func;
135: eps->which_ctx = data.which_ctx;
136: }
138: STGetMatMode(eps->OP,&matmode);
139: if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
140: /* Purify eigenvectors before reverting operator */
141: (*eps->ops->computevectors)(eps);
142: }
143: STPostSolve(eps->OP);
144: PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
146: if (!eps->reason) SETERRQ(((PetscObject)eps)->comm,1,"Internal error, solver returned without setting converged reason");
148: /* Map eigenvalues back to the original problem, necessary in some
149: * spectral transformations */
150: if (eps->ops->backtransform) {
151: (*eps->ops->backtransform)(eps);
152: }
154: /* Adjust left eigenvectors in generalized problems: y = B^T y */
155: if (eps->isgeneralized && eps->leftvecs) {
156: STGetOperators(eps->OP,PETSC_NULL,&B);
157: KSPCreate(((PetscObject)eps)->comm,&ksp);
158: KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
159: KSPSetFromOptions(ksp);
160: MatGetVecs(B,PETSC_NULL,&w);
161: for (i=0;i<eps->nconv;i++) {
162: VecCopy(eps->W[i],w);
163: KSPSolveTranspose(ksp,w,eps->W[i]);
164: }
165: KSPDestroy(&ksp);
166: VecDestroy(&w);
167: }
169: #if !defined(PETSC_USE_COMPLEX)
170: /* reorder conjugate eigenvalues (positive imaginary first) */
171: for (i=0; i<eps->nconv-1; i++) {
172: if (eps->eigi[i] != 0) {
173: if (eps->eigi[i] < 0) {
174: eps->eigi[i] = -eps->eigi[i];
175: eps->eigi[i+1] = -eps->eigi[i+1];
176: if (!eps->evecsavailable) {
177: /* the next correction only works with eigenvectors */
178: (*eps->ops->computevectors)(eps);
179: }
180: VecScale(eps->V[i+1],-1.0);
181: }
182: i++;
183: }
184: }
185: #endif
187: /* quick and dirty solution for FOLD: recompute eigenvalues as Rayleigh quotients */
188: PetscObjectTypeCompare((PetscObject)eps->OP,STFOLD,&isfold);
189: if (isfold) {
190: STGetOperators(eps->OP,&A,&B);
191: MatGetVecs(A,&w,PETSC_NULL);
192: if (!eps->evecsavailable) { (*eps->ops->computevectors)(eps); }
193: for (i=0;i<eps->nconv;i++) {
194: x = eps->V[i];
195: MatMult(A,x,w);
196: VecDot(w,x,&eps->eigr[i]);
197: if (eps->isgeneralized) {
198: MatMult(B,x,w);
199: VecDot(w,x,&dot);
200: eps->eigr[i] /= dot;
201: }
202: }
203: VecDestroy(&w);
204: }
206: /* In the case of Cayley transform, eigenvectors need to be B-normalized */
207: PetscObjectTypeCompare((PetscObject)eps->OP,STCAYLEY,&iscayley);
208: if (iscayley && eps->isgeneralized && eps->ishermitian) {
209: STGetOperators(eps->OP,PETSC_NULL,&B);
210: MatGetVecs(B,PETSC_NULL,&w);
211: if (!eps->evecsavailable) { (*eps->ops->computevectors)(eps); }
212: for (i=0;i<eps->nconv;i++) {
213: x = eps->V[i];
214: MatMult(B,x,w);
215: VecDot(w,x,&dot);
216: VecScale(x,1.0/PetscSqrtScalar(dot));
217: }
218: VecDestroy(&w);
219: }
221: /* sort eigenvalues according to eps->which parameter */
222: EPSSortEigenvalues(eps,eps->nconv,eps->eigr,eps->eigi,eps->perm);
224: if (!viewed) {
225: PetscOptionsGetString(((PetscObject)eps)->prefix,"-eps_view",filename,PETSC_MAX_PATH_LEN,&flg);
226: if (flg && !PetscPreLoadingOn) {
227: PetscViewerASCIIOpen(((PetscObject)eps)->comm,filename,&viewer);
228: EPSView(eps,viewer);
229: PetscViewerDestroy(&viewer);
230: }
231: }
233: flg = PETSC_FALSE;
234: PetscOptionsGetBool(((PetscObject)eps)->prefix,"-eps_plot_eigs",&flg,PETSC_NULL);
235: if (flg) {
236: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
237: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
238: PetscViewerDrawGetDraw(viewer,0,&draw);
239: PetscDrawSPCreate(draw,1,&drawsp);
240: for (i=0;i<eps->nconv;i++) {
241: #if defined(PETSC_USE_COMPLEX)
242: re = PetscRealPart(eps->eigr[i]);
243: im = PetscImaginaryPart(eps->eigi[i]);
244: #else
245: re = eps->eigr[i];
246: im = eps->eigi[i];
247: #endif
248: PetscDrawSPAddPoint(drawsp,&re,&im);
249: }
250: PetscDrawSPDraw(drawsp);
251: PetscDrawSPDestroy(&drawsp);
252: PetscViewerDestroy(&viewer);
253: }
255: /* Remove the initial subspaces */
256: eps->nini = 0;
257: eps->ninil = 0;
258: return(0);
259: }
263: /*@
264: EPSGetIterationNumber - Gets the current iteration number. If the
265: call to EPSSolve() is complete, then it returns the number of iterations
266: carried out by the solution method.
267:
268: Not Collective
270: Input Parameter:
271: . eps - the eigensolver context
273: Output Parameter:
274: . its - number of iterations
276: Level: intermediate
278: Note:
279: During the i-th iteration this call returns i-1. If EPSSolve() is
280: complete, then parameter "its" contains either the iteration number at
281: which convergence was successfully reached, or failure was detected.
282: Call EPSGetConvergedReason() to determine if the solver converged or
283: failed and why.
285: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
286: @*/
287: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
288: {
292: *its = eps->its;
293: return(0);
294: }
298: /*@
299: EPSGetOperationCounters - Gets the total number of operator applications,
300: inner product operations and linear iterations used by the ST object
301: during the last EPSSolve() call.
303: Not Collective
305: Input Parameter:
306: . eps - EPS context
308: Output Parameter:
309: + ops - number of operator applications
310: . dots - number of inner product operations
311: - lits - number of linear iterations
313: Notes:
314: When the eigensolver algorithm invokes STApply() then a linear system
315: must be solved (except in the case of standard eigenproblems and shift
316: transformation). The number of iterations required in this solve is
317: accumulated into a counter whose value is returned by this function.
319: These counters are reset to zero at each successive call to EPSSolve().
321: Level: intermediate
323: @*/
324: PetscErrorCode EPSGetOperationCounters(EPS eps,PetscInt* ops,PetscInt* dots,PetscInt* lits)
325: {
330: if (!eps->OP) { EPSGetST(eps,&eps->OP); }
331: STGetOperationCounters(eps->OP,ops,lits);
332: if (dots) {
333: if (!eps->ip) { EPSGetIP(eps,&eps->ip); }
334: IPGetOperationCounters(eps->ip,dots);
335: }
336: return(0);
337: }
341: /*@
342: EPSGetConverged - Gets the number of converged eigenpairs.
344: Not Collective
346: Input Parameter:
347: . eps - the eigensolver context
348:
349: Output Parameter:
350: . nconv - number of converged eigenpairs
352: Note:
353: This function should be called after EPSSolve() has finished.
355: Level: beginner
357: .seealso: EPSSetDimensions(), EPSSolve()
358: @*/
359: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
360: {
364: *nconv = eps->nconv;
365: return(0);
366: }
370: /*@C
371: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
372: stopped.
374: Not Collective
376: Input Parameter:
377: . eps - the eigensolver context
379: Output Parameter:
380: . reason - negative value indicates diverged, positive value converged
382: Possible values for reason:
383: + EPS_CONVERGED_TOL - converged up to tolerance
384: . EPS_DIVERGED_ITS - required more than its to reach convergence
385: - EPS_DIVERGED_BREAKDOWN - generic breakdown in method
387: Note:
388: Can only be called after the call to EPSSolve() is complete.
390: Level: intermediate
392: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
393: @*/
394: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
395: {
399: *reason = eps->reason;
400: return(0);
401: }
405: /*@
406: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
407: subspace.
409: Not Collective, but vectors are shared by all processors that share the EPS
411: Input Parameter:
412: . eps - the eigensolver context
413:
414: Output Parameter:
415: . v - an array of vectors
417: Notes:
418: This function should be called after EPSSolve() has finished.
420: The user should provide in v an array of nconv vectors, where nconv is
421: the value returned by EPSGetConverged().
423: The first k vectors returned in v span an invariant subspace associated
424: with the first k computed eigenvalues (note that this is not true if the
425: k-th eigenvalue is complex and matrix A is real; in this case the first
426: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
427: in X for all x in X (a similar definition applies for generalized
428: eigenproblems).
430: Level: intermediate
432: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspaceLeft()
433: @*/
434: PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec *v)
435: {
437: PetscInt i;
443: if (!eps->V) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
444: if (!eps->ishermitian && eps->evecsavailable) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector,EPSComputeRelativeError or EPSComputeResidualNorm");
445: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
446: for (i=0;i<eps->nconv;i++) {
447: VecPointwiseDivide(v[i],eps->V[i],eps->D);
448: VecNormalize(v[i],PETSC_NULL);
449: }
450: }
451: else {
452: for (i=0;i<eps->nconv;i++) {
453: VecCopy(eps->V[i],v[i]);
454: }
455: }
456: return(0);
457: }
461: /*@
462: EPSGetInvariantSubspaceLeft - Gets an orthonormal basis of the computed left
463: invariant subspace (only available in two-sided eigensolvers).
465: Not Collective, but vectors are shared by all processors that share the EPS
467: Input Parameter:
468: . eps - the eigensolver context
469:
470: Output Parameter:
471: . v - an array of vectors
473: Notes:
474: This function should be called after EPSSolve() has finished.
476: The user should provide in v an array of nconv vectors, where nconv is
477: the value returned by EPSGetConverged().
479: The first k vectors returned in v span a left invariant subspace associated
480: with the first k computed eigenvalues (note that this is not true if the
481: k-th eigenvalue is complex and matrix A is real; in this case the first
482: k+1 vectors should be used). A left invariant subspace Y of A satisfies y'A
483: in Y for all y in Y (a similar definition applies for generalized
484: eigenproblems).
486: Level: intermediate
488: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
489: @*/
490: PetscErrorCode EPSGetInvariantSubspaceLeft(EPS eps,Vec *v)
491: {
493: PetscInt i;
499: if (!eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with EPSSetLeftVectorsWanted");
500: if (!eps->W) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
501: if (!eps->ishermitian && eps->evecsavailable) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSGetInvariantSubspaceLeft must be called before EPSGetEigenpairLeft,EPSComputeRelativeErrorLeft or EPSComputeResidualNormLeft");
502: for (i=0;i<eps->nconv;i++) {
503: VecCopy(eps->W[i],v[i]);
504: }
505: return(0);
506: }
510: /*@
511: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
512: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
514: Logically Collective on EPS
516: Input Parameters:
517: + eps - eigensolver context
518: - i - index of the solution
520: Output Parameters:
521: + eigr - real part of eigenvalue
522: . eigi - imaginary part of eigenvalue
523: . Vr - real part of eigenvector
524: - Vi - imaginary part of eigenvector
526: Notes:
527: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
528: configured with complex scalars the eigenvalue is stored
529: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
530: set to zero).
532: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
533: Eigenpairs are indexed according to the ordering criterion established
534: with EPSSetWhichEigenpairs().
536: The 2-norm of the eigenvector is one unless the problem is generalized
537: Hermitian. In this case the eigenvector is normalized with respect to the
538: norm defined by the B matrix.
540: Level: beginner
542: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSGetEigenvectorLeft(), EPSSolve(),
543: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
544: @*/
545: PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
546: {
552: if (!eps->eigr || !eps->eigi || !eps->V) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
553: if (i<0 || i>=eps->nconv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
554: EPSGetEigenvalue(eps,i,eigr,eigi);
555: if (Vr || Vi) { EPSGetEigenvector(eps,i,Vr,Vi); }
556: return(0);
557: }
561: /*@
562: EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve().
564: Not Collective
566: Input Parameters:
567: + eps - eigensolver context
568: - i - index of the solution
570: Output Parameters:
571: + eigr - real part of eigenvalue
572: - eigi - imaginary part of eigenvalue
574: Notes:
575: If the eigenvalue is real, then eigi is set to zero. If PETSc is
576: configured with complex scalars the eigenvalue is stored
577: directly in eigr (eigi is set to zero).
579: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
580: Eigenpairs are indexed according to the ordering criterion established
581: with EPSSetWhichEigenpairs().
583: Level: beginner
585: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
586: EPSGetEigenpair()
587: @*/
588: PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi)
589: {
590: PetscInt k;
594: if (!eps->eigr || !eps->eigi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
595: if (i<0 || i>=eps->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
596: if (!eps->perm) k = i;
597: else k = eps->perm[i];
598: #if defined(PETSC_USE_COMPLEX)
599: if (eigr) *eigr = eps->eigr[k];
600: if (eigi) *eigi = 0;
601: #else
602: if (eigr) *eigr = eps->eigr[k];
603: if (eigi) *eigi = eps->eigi[k];
604: #endif
605: return(0);
606: }
610: /*@
611: EPSGetEigenvector - Gets the i-th right eigenvector as computed by EPSSolve().
613: Logically Collective on EPS
615: Input Parameters:
616: + eps - eigensolver context
617: - i - index of the solution
619: Output Parameters:
620: + Vr - real part of eigenvector
621: - Vi - imaginary part of eigenvector
623: Notes:
624: If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
625: configured with complex scalars the eigenvector is stored
626: directly in Vr (Vi is set to zero).
628: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
629: Eigenpairs are indexed according to the ordering criterion established
630: with EPSSetWhichEigenpairs().
632: The 2-norm of the eigenvector is one unless the problem is generalized
633: Hermitian. In this case the eigenvector is normalized with respect to the
634: norm defined by the B matrix.
636: Level: beginner
638: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
639: EPSGetEigenpair(), EPSGetEigenvectorLeft()
640: @*/
641: PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi)
642: {
644: PetscInt k;
652: if (!eps->V) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
653: if (i<0 || i>=eps->nconv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
654: if (!eps->evecsavailable) {
655: (*eps->ops->computevectors)(eps);
656: }
657: if (!eps->perm) k = i;
658: else k = eps->perm[i];
659: #if defined(PETSC_USE_COMPLEX)
660: VecCopy(eps->V[k],Vr);
661: if (Vi) { VecSet(Vi,0.0); }
662: #else
663: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
664: VecCopy(eps->V[k],Vr);
665: if (Vi) { VecCopy(eps->V[k+1],Vi); }
666: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
667: VecCopy(eps->V[k-1],Vr);
668: if (Vi) {
669: VecCopy(eps->V[k],Vi);
670: VecScale(Vi,-1.0);
671: }
672: } else { /* real eigenvalue */
673: VecCopy(eps->V[k],Vr);
674: if (Vi) { VecSet(Vi,0.0); }
675: }
676: #endif
677: return(0);
678: }
682: /*@
683: EPSGetEigenvectorLeft - Gets the i-th left eigenvector as computed by EPSSolve()
684: (only available in two-sided eigensolvers).
686: Logically Collective on EPS
688: Input Parameters:
689: + eps - eigensolver context
690: - i - index of the solution
692: Output Parameters:
693: + Wr - real part of eigenvector
694: - Wi - imaginary part of eigenvector
696: Notes:
697: If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is
698: configured with complex scalars the eigenvector is stored
699: directly in Wr (Wi is set to zero).
701: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
702: Eigenpairs are indexed according to the ordering criterion established
703: with EPSSetWhichEigenpairs().
705: Level: beginner
707: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(),
708: EPSGetEigenpair(), EPSGetEigenvector()
709: @*/
710: PetscErrorCode EPSGetEigenvectorLeft(EPS eps,PetscInt i,Vec Wr,Vec Wi)
711: {
713: PetscInt k;
721: if (!eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with EPSSetLeftVectorsWanted");
722: if (!eps->W) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
723: if (i<0 || i>=eps->nconv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
724: if (!eps->evecsavailable) {
725: (*eps->ops->computevectors)(eps);
726: }
727: if (!eps->perm) k = i;
728: else k = eps->perm[i];
729: #if defined(PETSC_USE_COMPLEX)
730: VecCopy(eps->W[k],Wr);
731: if (Wi) { VecSet(Wi,0.0); }
732: #else
733: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
734: VecCopy(eps->W[k],Wr);
735: if (Wi) { VecCopy(eps->W[k+1],Wi); }
736: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
737: VecCopy(eps->W[k-1],Wr);
738: if (Wi) {
739: VecCopy(eps->W[k],Wi);
740: VecScale(Wi,-1.0);
741: }
742: } else { /* real eigenvalue */
743: VecCopy(eps->W[k],Wr);
744: if (Wi) { VecSet(Wi,0.0); }
745: }
746: #endif
747: return(0);
748: }
752: /*@
753: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
754: computed eigenpair.
756: Not Collective
758: Input Parameter:
759: + eps - eigensolver context
760: - i - index of eigenpair
762: Output Parameter:
763: . errest - the error estimate
765: Notes:
766: This is the error estimate used internally by the eigensolver. The actual
767: error bound can be computed with EPSComputeRelativeError(). See also the users
768: manual for details.
770: Level: advanced
772: .seealso: EPSComputeRelativeError()
773: @*/
774: PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)
775: {
779: if (!eps->eigr || !eps->eigi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
780: if (i<0 || i>=eps->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
781: if (eps->perm) i = eps->perm[i];
782: if (errest) *errest = eps->errest[i];
783: return(0);
784: }
788: /*@
789: EPSGetErrorEstimateLeft - Returns the left error estimate associated to the i-th
790: computed eigenpair (only available in two-sided eigensolvers).
792: Not Collective
794: Input Parameter:
795: + eps - eigensolver context
796: - i - index of eigenpair
798: Output Parameter:
799: . errest - the left error estimate
801: Notes:
802: This is the error estimate used internally by the eigensolver. The actual
803: error bound can be computed with EPSComputeRelativeErrorLeft(). See also the users
804: manual for details.
806: Level: advanced
808: .seealso: EPSComputeRelativeErrorLeft()
809: @*/
810: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps,PetscInt i,PetscReal *errest)
811: {
815: if (!eps->eigr || !eps->eigi) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"EPSSolve must be called first");
816: if (!eps->leftvecs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with EPSSetLeftVectorsWanted");
817: if (i<0 || i>=eps->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
818: if (eps->perm) i = eps->perm[i];
819: if (errest) *errest = eps->errest_left[i];
820: return(0);
821: }
825: /*
826: EPSComputeResidualNorm_Private - Computes the norm of the residual vector
827: associated with an eigenpair.
828: */
829: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *norm)
830: {
832: Vec u,w;
833: Mat A,B;
834: #if !defined(PETSC_USE_COMPLEX)
835: Vec v;
836: PetscReal ni,nr;
837: #endif
838:
840: STGetOperators(eps->OP,&A,&B);
841: VecDuplicate(eps->V[0],&u);
842: VecDuplicate(eps->V[0],&w);
843:
844: #if !defined(PETSC_USE_COMPLEX)
845: if (ki == 0 ||
846: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
847: #endif
848: MatMult(A,xr,u); /* u=A*x */
849: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
850: if (eps->isgeneralized) { MatMult(B,xr,w); }
851: else { VecCopy(xr,w); } /* w=B*x */
852: VecAXPY(u,-kr,w); /* u=A*x-k*B*x */
853: }
854: VecNorm(u,NORM_2,norm);
855: #if !defined(PETSC_USE_COMPLEX)
856: } else {
857: VecDuplicate(eps->V[0],&v);
858: MatMult(A,xr,u); /* u=A*xr */
859: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
860: if (eps->isgeneralized) { MatMult(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) { MatMult(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: }
867: VecNorm(u,NORM_2,&nr);
868: MatMult(A,xi,u); /* u=A*xi */
869: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
870: VecAXPY(u,-kr,w); /* u=A*xi-kr*B*xi */
871: VecAXPY(u,-ki,v); /* u=A*xi-kr*B*xi-ki*B*xr */
872: }
873: VecNorm(u,NORM_2,&ni);
874: *norm = SlepcAbsEigenvalue(nr,ni);
875: VecDestroy(&v);
876: }
877: #endif
879: VecDestroy(&w);
880: VecDestroy(&u);
881: return(0);
882: }
886: /*@
887: EPSComputeResidualNorm - Computes the norm of the residual vector associated with
888: the i-th computed eigenpair.
890: Collective on EPS
892: Input Parameter:
893: . eps - the eigensolver context
894: . i - the solution index
896: Output Parameter:
897: . norm - the residual norm, computed as ||Ax-kBx||_2 where k is the
898: eigenvalue and x is the eigenvector.
899: If k=0 then the residual norm is computed as ||Ax||_2.
901: Notes:
902: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
903: Eigenpairs are indexed according to the ordering criterion established
904: with EPSSetWhichEigenpairs().
906: Level: beginner
908: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
909: @*/
910: PetscErrorCode EPSComputeResidualNorm(EPS eps,PetscInt i,PetscReal *norm)
911: {
913: Vec xr,xi;
914: PetscScalar kr,ki;
915:
920: VecDuplicate(eps->V[0],&xr);
921: VecDuplicate(eps->V[0],&xi);
922: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
923: EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,norm);
924: VecDestroy(&xr);
925: VecDestroy(&xi);
926: return(0);
927: }
931: /*@
932: EPSComputeResidualNormLeft - Computes the norm of the residual vector associated with
933: the i-th computed left eigenvector (only available in two-sided eigensolvers).
935: Collective on EPS
937: Input Parameter:
938: . eps - the eigensolver context
939: . i - the solution index
941: Output Parameter:
942: . norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the
943: eigenvalue and y is the left eigenvector.
944: If k=0 then the residual norm is computed as ||y'A||_2.
946: Notes:
947: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
948: Eigenpairs are indexed according to the ordering criterion established
949: with EPSSetWhichEigenpairs().
951: Level: beginner
953: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
954: @*/
955: PetscErrorCode EPSComputeResidualNormLeft(EPS eps,PetscInt i,PetscReal *norm)
956: {
958: Vec u,v,w,xr,xi;
959: Mat A,B;
960: PetscScalar kr,ki;
961: #if !defined(PETSC_USE_COMPLEX)
962: PetscReal ni,nr;
963: #endif
964:
969: if (!eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with EPSSetLeftVectorsWanted");
970: STGetOperators(eps->OP,&A,&B);
971: VecDuplicate(eps->W[0],&u);
972: VecDuplicate(eps->W[0],&v);
973: VecDuplicate(eps->W[0],&w);
974: VecDuplicate(eps->W[0],&xr);
975: VecDuplicate(eps->W[0],&xi);
976: EPSGetEigenvalue(eps,i,&kr,&ki);
977: EPSGetEigenvectorLeft(eps,i,xr,xi);
979: #if !defined(PETSC_USE_COMPLEX)
980: if (ki == 0 ||
981: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
982: #endif
983: MatMultTranspose(A,xr,u); /* u=A'*x */
984: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
985: if (eps->isgeneralized) { MatMultTranspose(B,xr,w); }
986: else { VecCopy(xr,w); } /* w=B'*x */
987: VecAXPY(u,-kr,w); /* u=A'*x-k*B'*x */
988: }
989: VecNorm(u,NORM_2,norm);
990: #if !defined(PETSC_USE_COMPLEX)
991: } else {
992: MatMultTranspose(A,xr,u); /* u=A'*xr */
993: if (eps->isgeneralized) { MatMultTranspose(B,xr,v); }
994: else { VecCopy(xr,v); } /* v=B'*xr */
995: VecAXPY(u,-kr,v); /* u=A'*xr-kr*B'*xr */
996: if (eps->isgeneralized) { MatMultTranspose(B,xi,w); }
997: else { VecCopy(xi,w); } /* w=B'*xi */
998: VecAXPY(u,ki,w); /* u=A'*xr-kr*B'*xr+ki*B'*xi */
999: VecNorm(u,NORM_2,&nr);
1000: MatMultTranspose(A,xi,u); /* u=A'*xi */
1001: VecAXPY(u,-kr,w); /* u=A'*xi-kr*B'*xi */
1002: VecAXPY(u,-ki,v); /* u=A'*xi-kr*B'*xi-ki*B'*xr */
1003: VecNorm(u,NORM_2,&ni);
1004: *norm = SlepcAbsEigenvalue(nr,ni);
1005: }
1006: #endif
1008: VecDestroy(&w);
1009: VecDestroy(&v);
1010: VecDestroy(&u);
1011: VecDestroy(&xr);
1012: VecDestroy(&xi);
1013: return(0);
1014: }
1018: /*
1019: EPSComputeRelativeError_Private - Computes the relative error bound
1020: associated with an eigenpair.
1021: */
1022: PetscErrorCode EPSComputeRelativeError_Private(EPS eps,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *error)
1023: {
1025: PetscReal norm,er;
1026: #if !defined(PETSC_USE_COMPLEX)
1027: PetscReal ei;
1028: #endif
1029:
1031: EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,&norm);
1033: #if !defined(PETSC_USE_COMPLEX)
1034: if (ki == 0) {
1035: #endif
1036: VecNorm(xr,NORM_2,&er);
1037: #if !defined(PETSC_USE_COMPLEX)
1038: } else {
1039: VecNorm(xr,NORM_2,&er);
1040: VecNorm(xi,NORM_2,&ei);
1041: er = SlepcAbsEigenvalue(er,ei);
1042: }
1043: #endif
1044: (*eps->conv_func)(eps,kr,ki,norm/er,error,eps->conv_ctx);
1045: return(0);
1046: }
1050: /*@
1051: EPSComputeRelativeError - Computes the relative error bound associated
1052: with the i-th computed eigenpair.
1054: Collective on EPS
1056: Input Parameter:
1057: . eps - the eigensolver context
1058: . i - the solution index
1060: Output Parameter:
1061: . error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where
1062: k is the eigenvalue and x is the eigenvector.
1063: If k=0 the relative error is computed as ||Ax||_2/||x||_2.
1065: Level: beginner
1067: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
1068: @*/
1069: PetscErrorCode EPSComputeRelativeError(EPS eps,PetscInt i,PetscReal *error)
1070: {
1072: Vec xr,xi;
1073: PetscScalar kr,ki;
1074:
1079: VecDuplicate(eps->V[0],&xr);
1080: VecDuplicate(eps->V[0],&xi);
1081: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
1082: EPSComputeRelativeError_Private(eps,kr,ki,xr,xi,error);
1083: VecDestroy(&xr);
1084: VecDestroy(&xi);
1085: return(0);
1086: }
1090: /*@
1091: EPSComputeRelativeErrorLeft - Computes the relative error bound associated
1092: with the i-th computed eigenvalue and left eigenvector (only available in
1093: two-sided eigensolvers).
1095: Collective on EPS
1097: Input Parameter:
1098: . eps - the eigensolver context
1099: . i - the solution index
1101: Output Parameter:
1102: . error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where
1103: k is the eigenvalue and y is the left eigenvector.
1104: If k=0 the relative error is computed as ||y'A||_2/||y||_2.
1106: Level: beginner
1108: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
1109: @*/
1110: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps,PetscInt i,PetscReal *error)
1111: {
1113: Vec xr,xi;
1114: PetscScalar kr,ki;
1115: PetscReal norm,er;
1116: #if !defined(PETSC_USE_COMPLEX)
1117: Vec u;
1118: PetscReal ei;
1119: #endif
1120:
1123: EPSComputeResidualNormLeft(eps,i,&norm);
1126: VecDuplicate(eps->W[0],&xr);
1127: VecDuplicate(eps->W[0],&xi);
1128: EPSGetEigenvalue(eps,i,&kr,&ki);
1129: EPSGetEigenvectorLeft(eps,i,xr,xi);
1131: #if !defined(PETSC_USE_COMPLEX)
1132: if (ki == 0 ||
1133: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1134: #endif
1135: VecNorm(xr,NORM_2,&er);
1136: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1137: *error = norm / (PetscAbsScalar(kr) * er);
1138: } else {
1139: *error = norm / er;
1140: }
1141: #if !defined(PETSC_USE_COMPLEX)
1142: } else {
1143: VecDuplicate(xi,&u);
1144: VecCopy(xi,u);
1145: VecAXPBY(u,kr,-ki,xr);
1146: VecNorm(u,NORM_2,&er);
1147: VecAXPBY(xi,kr,ki,xr);
1148: VecNorm(xi,NORM_2,&ei);
1149: VecDestroy(&u);
1150: *error = norm / SlepcAbsEigenvalue(er,ei);
1151: }
1152: #endif
1153:
1154: VecDestroy(&xr);
1155: VecDestroy(&xi);
1156: return(0);
1157: }
1159: #define SWAP(a,b,t) {t=a;a=b;b=t;}
1163: /*@
1164: EPSSortEigenvalues - Sorts a list of eigenvalues according to the criterion
1165: specified via EPSSetWhichEigenpairs().
1167: Not Collective
1169: Input Parameters:
1170: + eps - the eigensolver context
1171: . n - number of eigenvalues in the list
1172: . eigr - pointer to the array containing the eigenvalues
1173: - eigi - imaginary part of the eigenvalues (only when using real numbers)
1175: Output Parameter:
1176: . perm - resulting permutation
1178: Note:
1179: The result is a list of indices in the original eigenvalue array
1180: corresponding to the first nev eigenvalues sorted in the specified
1181: criterion.
1183: Level: developer
1185: .seealso: EPSSetWhichEigenpairs()
1186: @*/
1187: PetscErrorCode EPSSortEigenvalues(EPS eps,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
1188: {
1190: PetscScalar re,im;
1191: PetscInt i,j,result,tmp;
1198: for (i=0; i<n; i++) { perm[i] = i; }
1199: /* insertion sort */
1200: for (i=n-1; i>=0; i--) {
1201: re = eigr[perm[i]];
1202: im = eigi[perm[i]];
1203: j = i + 1;
1204: #if !defined(PETSC_USE_COMPLEX)
1205: if (im != 0) {
1206: /* complex eigenvalue */
1207: i--;
1208: im = eigi[perm[i]];
1209: }
1210: #endif
1211: while (j<n) {
1212: EPSCompareEigenvalues(eps,re,im,eigr[perm[j]],eigi[perm[j]],&result);
1213: if (result < 0) break;
1214: #if !defined(PETSC_USE_COMPLEX)
1215: /* keep together every complex conjugated eigenpair */
1216: if (im == 0) {
1217: if (eigi[perm[j]] == 0) {
1218: #endif
1219: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
1220: j++;
1221: #if !defined(PETSC_USE_COMPLEX)
1222: } else {
1223: tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
1224: j+=2;
1225: }
1226: } else {
1227: if (eigi[perm[j]] == 0) {
1228: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
1229: j++;
1230: } else {
1231: tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
1232: tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
1233: j+=2;
1234: }
1235: }
1236: #endif
1237: }
1238: }
1239: return(0);
1240: }
1244: /*@
1245: EPSCompareEigenvalues - Compares two (possibly complex) eigenvalues according
1246: to a certain criterion.
1248: Not Collective
1250: Input Parameters:
1251: + eps - the eigensolver context
1252: . ar - real part of the 1st eigenvalue
1253: . ai - imaginary part of the 1st eigenvalue
1254: . br - real part of the 2nd eigenvalue
1255: - bi - imaginary part of the 2nd eigenvalue
1257: Output Parameter:
1258: . res - result of comparison
1260: Notes:
1261: The returning parameter 'res' can be:
1262: + negative - if the 1st eigenvalue is preferred to the 2st one
1263: . zero - if both eigenvalues are equally preferred
1264: - positive - if the 2st eigenvalue is preferred to the 1st one
1266: The criterion of comparison is related to the 'which' parameter set with
1267: EPSSetWhichEigenpairs().
1269: Level: developer
1271: .seealso: EPSSortEigenvalues(), EPSSetWhichEigenpairs()
1272: @*/
1273: PetscErrorCode EPSCompareEigenvalues(EPS eps,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
1274: {
1280: if (!eps->which_func) SETERRQ(PETSC_COMM_SELF,1,"Undefined eigenvalue comparison function");
1281: (*eps->which_func)(ar,ai,br,bi,result,eps->which_ctx);
1282: return(0);
1283: }
1287: /*@
1288: EPSGetStartVector - Gets a suitable vector to be used as the starting vector
1289: for the recurrence that builds the right subspace.
1291: Collective on EPS and Vec
1293: Input Parameters:
1294: + eps - the eigensolver context
1295: - i - iteration number
1297: Output Parameters:
1298: + vec - the start vector
1299: - breakdown - flag indicating that a breakdown has occurred
1301: Notes:
1302: The start vector is computed from another vector: for the first step (i=0),
1303: the first initial vector is used (see EPSSetInitialSpace()); otherwise a random
1304: vector is created. Then this vector is forced to be in the range of OP (only
1305: for generalized definite problems) and orthonormalized with respect to all
1306: V-vectors up to i-1.
1308: The flag breakdown is set to true if either i=0 and the vector belongs to the
1309: deflation space, or i>0 and the vector is linearly dependent with respect
1310: to the V-vectors.
1312: The caller must pass a vector already allocated with dimensions conforming
1313: to the initial vector. This vector is overwritten.
1315: Level: developer
1317: .seealso: EPSSetInitialSpace()
1318: @*/
1319: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,Vec vec,PetscBool *breakdown)
1320: {
1322: PetscReal norm;
1323: PetscBool lindep;
1324: Vec w;
1325:
1332: VecDuplicate(eps->V[0],&w);
1334: /* For the first step, use the first initial vector, otherwise a random one */
1335: if (i==0 && eps->nini>0) {
1336: VecCopy(eps->V[0],w);
1337: } else {
1338: SlepcVecSetRandom(w,eps->rand);
1339: }
1341: /* Force the vector to be in the range of OP for definite generalized problems */
1342: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
1343: STApply(eps->OP,w,vec);
1344: } else {
1345: VecCopy(w,vec);
1346: }
1348: /* Orthonormalize the vector with respect to previous vectors */
1349: IPOrthogonalize(eps->ip,eps->nds,eps->defl,i,PETSC_NULL,eps->V,vec,PETSC_NULL,&norm,&lindep);
1350: if (breakdown) *breakdown = lindep;
1351: else if (lindep || norm == 0.0) {
1352: if (i==0) SETERRQ(((PetscObject)eps)->comm,1,"Initial vector is zero or belongs to the deflation space");
1353: else SETERRQ(((PetscObject)eps)->comm,1,"Unable to generate more start vectors");
1354: }
1355: VecScale(vec,1.0/norm);
1357: VecDestroy(&w);
1358: return(0);
1359: }
1363: /*@
1364: EPSGetStartVectorLeft - Gets a suitable vector to be used as the starting vector
1365: in the recurrence that builds the left subspace (in methods that work with two
1366: subspaces).
1368: Collective on EPS and Vec
1370: Input Parameters:
1371: + eps - the eigensolver context
1372: - i - iteration number
1374: Output Parameter:
1375: + vec - the start vector
1376: - breakdown - flag indicating that a breakdown has occurred
1378: Notes:
1379: The start vector is computed from another vector: for the first step (i=0),
1380: the first left initial vector is used (see EPSSetInitialSpaceLeft()); otherwise
1381: a random vector is created. Then this vector is forced to be in the range
1382: of OP' and orthonormalized with respect to all W-vectors up to i-1.
1384: The flag breakdown is set to true if i>0 and the vector is linearly dependent
1385: with respect to the W-vectors.
1387: The caller must pass a vector already allocated with dimensions conforming
1388: to the left initial vector. This vector is overwritten.
1390: Level: developer
1392: .seealso: EPSSetInitialSpaceLeft()
1394: @*/
1395: PetscErrorCode EPSGetStartVectorLeft(EPS eps,PetscInt i,Vec vec,PetscBool *breakdown)
1396: {
1398: PetscReal norm;
1399: PetscBool lindep;
1400: Vec w;
1401:
1408: VecDuplicate(eps->W[0],&w);
1410: /* For the first step, use the first initial left vector, otherwise a random one */
1411: if (i==0 && eps->ninil>0) {
1412: VecCopy(eps->W[0],w);
1413: } else {
1414: SlepcVecSetRandom(w,eps->rand);
1415: }
1417: /* Force the vector to be in the range of OP' */
1418: STApplyTranspose(eps->OP,w,vec);
1420: /* Orthonormalize the vector with respect to previous vectors */
1421: IPOrthogonalize(eps->ip,0,PETSC_NULL,i,PETSC_NULL,eps->W,vec,PETSC_NULL,&norm,&lindep);
1422: if (breakdown) *breakdown = lindep;
1423: else if (lindep || norm == 0.0) {
1424: if (i==0) SETERRQ(((PetscObject)eps)->comm,1,"Left initial vector is zero");
1425: else SETERRQ(((PetscObject)eps)->comm,1,"Unable to generate more left start vectors");
1426: }
1427: VecScale(vec,1/norm);
1429: VecDestroy(&w);
1430: return(0);
1431: }