1: /*
2: EPS routines related to the solution process.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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*/
25: #include <petscdraw.h>
29: PETSC_STATIC_INLINE PetscErrorCode EPSComputeVectors(EPS eps) 30: {
34: EPSCheckSolved(eps,1);
35: switch (eps->state) {
36: case EPS_STATE_SOLVED:
37: if (eps->ops->computevectors) {
38: (*eps->ops->computevectors)(eps);
39: }
40: break;
41: default: 42: break;
43: }
44: eps->state = EPS_STATE_EIGENVECTORS;
45: return(0);
46: }
50: /*@
51: EPSSolve - Solves the eigensystem.
53: Collective on EPS 55: Input Parameter:
56: . eps - eigensolver context obtained from EPSCreate()
58: Options Database Keys:
59: + -eps_view - print information about the solver used
60: . -eps_view_mat0 binary - save the first matrix (A) to the default binary viewer
61: . -eps_view_mat1 binary - save the second matrix (B) to the default binary viewer
62: - -eps_plot_eigs - plot computed eigenvalues
64: Level: beginner
66: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
67: @*/
68: PetscErrorCode EPSSolve(EPS eps) 69: {
70: PetscErrorCode ierr;
71: PetscInt i,nmat;
72: PetscReal re,im;
73: PetscScalar dot;
74: PetscBool flg,iscayley;
75: PetscViewer viewer;
76: PetscViewerFormat format;
77: PetscDraw draw;
78: PetscDrawSP drawsp;
79: STMatMode matmode;
80: Mat A,B;
81: Vec w,x;
85: PetscLogEventBegin(EPS_Solve,eps,0,0,0);
87: /* call setup */
88: EPSSetUp(eps);
89: eps->nconv = 0;
90: eps->its = 0;
91: for (i=0;i<eps->ncv;i++) {
92: eps->eigr[i] = 0.0;
93: eps->eigi[i] = 0.0;
94: eps->errest[i] = 0.0;
95: }
96: EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->ncv);
98: /* call solver */
99: (*eps->ops->solve)(eps);
100: eps->state = EPS_STATE_SOLVED;
102: STGetMatMode(eps->st,&matmode);
103: if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
104: /* Purify eigenvectors before reverting operator */
105: EPSComputeVectors(eps);
106: }
107: STPostSolve(eps->st);
109: if (!eps->reason) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
111: /* Map eigenvalues back to the original problem, necessary in some
112: * spectral transformations */
113: if (eps->ops->backtransform) {
114: (*eps->ops->backtransform)(eps);
115: }
117: #if !defined(PETSC_USE_COMPLEX)
118: /* reorder conjugate eigenvalues (positive imaginary first) */
119: for (i=0; i<eps->nconv-1; i++) {
120: if (eps->eigi[i] != 0) {
121: if (eps->eigi[i] < 0) {
122: eps->eigi[i] = -eps->eigi[i];
123: eps->eigi[i+1] = -eps->eigi[i+1];
124: /* the next correction only works with eigenvectors */
125: EPSComputeVectors(eps);
126: BVScaleColumn(eps->V,i+1,-1.0);
127: }
128: i++;
129: }
130: }
131: #endif
133: STGetNumMatrices(eps->st,&nmat);
134: STGetOperators(eps->st,0,&A);
135: if (nmat>1) { STGetOperators(eps->st,1,&B); }
137: /* In the case of Cayley transform, eigenvectors need to be B-normalized */
138: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
139: if (iscayley && eps->isgeneralized && eps->ishermitian) {
140: MatGetVecs(B,NULL,&w);
141: EPSComputeVectors(eps);
142: for (i=0;i<eps->nconv;i++) {
143: BVGetColumn(eps->V,i,&x);
144: MatMult(B,x,w);
145: VecDot(w,x,&dot);
146: VecScale(x,1.0/PetscSqrtScalar(dot));
147: BVRestoreColumn(eps->V,i,&x);
148: }
149: VecDestroy(&w);
150: }
152: /* sort eigenvalues according to eps->which parameter */
153: SlepcSortEigenvalues(eps->sc,eps->nconv,eps->eigr,eps->eigi,eps->perm);
155: PetscLogEventEnd(EPS_Solve,eps,0,0,0);
157: /* various viewers */
158: MatViewFromOptions(A,((PetscObject)eps)->prefix,"-eps_view_mat0");
159: if (nmat>1) { MatViewFromOptions(B,((PetscObject)eps)->prefix,"-eps_view_mat1"); }
161: PetscOptionsGetViewer(PetscObjectComm((PetscObject)eps),((PetscObject)eps)->prefix,"-eps_view",&viewer,&format,&flg);
162: if (flg && !PetscPreLoadingOn) {
163: PetscViewerPushFormat(viewer,format);
164: EPSView(eps,viewer);
165: PetscViewerPopFormat(viewer);
166: PetscViewerDestroy(&viewer);
167: }
169: flg = PETSC_FALSE;
170: PetscOptionsGetBool(((PetscObject)eps)->prefix,"-eps_plot_eigs",&flg,NULL);
171: if (flg) {
172: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
173: PetscViewerDrawGetDraw(viewer,0,&draw);
174: PetscDrawSPCreate(draw,1,&drawsp);
175: for (i=0;i<eps->nconv;i++) {
176: #if defined(PETSC_USE_COMPLEX)
177: re = PetscRealPart(eps->eigr[i]);
178: im = PetscImaginaryPart(eps->eigi[i]);
179: #else
180: re = eps->eigr[i];
181: im = eps->eigi[i];
182: #endif
183: PetscDrawSPAddPoint(drawsp,&re,&im);
184: }
185: PetscDrawSPDraw(drawsp,PETSC_TRUE);
186: PetscDrawSPDestroy(&drawsp);
187: PetscViewerDestroy(&viewer);
188: }
190: /* Remove deflation and initial subspaces */
191: eps->nds = 0;
192: eps->nini = 0;
193: return(0);
194: }
198: /*@
199: EPSGetIterationNumber - Gets the current iteration number. If the
200: call to EPSSolve() is complete, then it returns the number of iterations
201: carried out by the solution method.
203: Not Collective
205: Input Parameter:
206: . eps - the eigensolver context
208: Output Parameter:
209: . its - number of iterations
211: Level: intermediate
213: Note:
214: During the i-th iteration this call returns i-1. If EPSSolve() is
215: complete, then parameter "its" contains either the iteration number at
216: which convergence was successfully reached, or failure was detected.
217: Call EPSGetConvergedReason() to determine if the solver converged or
218: failed and why.
220: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
221: @*/
222: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)223: {
227: *its = eps->its;
228: return(0);
229: }
233: /*@
234: EPSGetConverged - Gets the number of converged eigenpairs.
236: Not Collective
238: Input Parameter:
239: . eps - the eigensolver context
241: Output Parameter:
242: . nconv - number of converged eigenpairs
244: Note:
245: This function should be called after EPSSolve() has finished.
247: Level: beginner
249: .seealso: EPSSetDimensions(), EPSSolve()
250: @*/
251: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)252: {
256: EPSCheckSolved(eps,1);
257: *nconv = eps->nconv;
258: return(0);
259: }
263: /*@C
264: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
265: stopped.
267: Not Collective
269: Input Parameter:
270: . eps - the eigensolver context
272: Output Parameter:
273: . reason - negative value indicates diverged, positive value converged
275: Possible values for reason:
276: + EPS_CONVERGED_TOL - converged up to tolerance
277: . EPS_DIVERGED_ITS - required more than its to reach convergence
278: - EPS_DIVERGED_BREAKDOWN - generic breakdown in method
280: Note:
281: Can only be called after the call to EPSSolve() is complete.
283: Level: intermediate
285: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason286: @*/
287: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)288: {
292: EPSCheckSolved(eps,1);
293: *reason = eps->reason;
294: return(0);
295: }
299: /*@
300: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
301: subspace.
303: Not Collective, but vectors are shared by all processors that share the EPS305: Input Parameter:
306: . eps - the eigensolver context
308: Output Parameter:
309: . v - an array of vectors
311: Notes:
312: This function should be called after EPSSolve() has finished.
314: The user should provide in v an array of nconv vectors, where nconv is
315: the value returned by EPSGetConverged().
317: The first k vectors returned in v span an invariant subspace associated
318: with the first k computed eigenvalues (note that this is not true if the
319: k-th eigenvalue is complex and matrix A is real; in this case the first
320: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
321: in X for all x in X (a similar definition applies for generalized
322: eigenproblems).
324: Level: intermediate
326: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve()
327: @*/
328: PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec *v)329: {
331: PetscInt i;
337: EPSCheckSolved(eps,1);
338: if (!eps->ishermitian && eps->state==EPS_STATE_EIGENVECTORS) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector,EPSComputeRelativeError or EPSComputeResidualNorm");
339: for (i=0;i<eps->nconv;i++) {
340: BVCopyVec(eps->V,i,v[i]);
341: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
342: VecPointwiseDivide(v[i],v[i],eps->D);
343: VecNormalize(v[i],NULL);
344: }
345: }
346: return(0);
347: }
351: /*@
352: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
353: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
355: Logically Collective on EPS357: Input Parameters:
358: + eps - eigensolver context
359: - i - index of the solution
361: Output Parameters:
362: + eigr - real part of eigenvalue
363: . eigi - imaginary part of eigenvalue
364: . Vr - real part of eigenvector
365: - Vi - imaginary part of eigenvector
367: Notes:
368: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
369: configured with complex scalars the eigenvalue is stored
370: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
371: set to zero).
373: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
374: Eigenpairs are indexed according to the ordering criterion established
375: with EPSSetWhichEigenpairs().
377: The 2-norm of the eigenvector is one unless the problem is generalized
378: Hermitian. In this case the eigenvector is normalized with respect to the
379: norm defined by the B matrix.
381: Level: beginner
383: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSSolve(),
384: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
385: @*/
386: PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)387: {
393: EPSCheckSolved(eps,1);
394: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
395: EPSGetEigenvalue(eps,i,eigr,eigi);
396: if (Vr || Vi) { EPSGetEigenvector(eps,i,Vr,Vi); }
397: return(0);
398: }
402: /*@
403: EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve().
405: Not Collective
407: Input Parameters:
408: + eps - eigensolver context
409: - i - index of the solution
411: Output Parameters:
412: + eigr - real part of eigenvalue
413: - eigi - imaginary part of eigenvalue
415: Notes:
416: If the eigenvalue is real, then eigi is set to zero. If PETSc is
417: configured with complex scalars the eigenvalue is stored
418: directly in eigr (eigi is set to zero).
420: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
421: Eigenpairs are indexed according to the ordering criterion established
422: with EPSSetWhichEigenpairs().
424: Level: beginner
426: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
427: @*/
428: PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi)429: {
430: PetscInt k;
434: EPSCheckSolved(eps,1);
435: if (i<0 || i>=eps->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
436: if (!eps->perm) k = i;
437: else k = eps->perm[i];
438: #if defined(PETSC_USE_COMPLEX)
439: if (eigr) *eigr = eps->eigr[k];
440: if (eigi) *eigi = 0;
441: #else
442: if (eigr) *eigr = eps->eigr[k];
443: if (eigi) *eigi = eps->eigi[k];
444: #endif
445: return(0);
446: }
450: /*@
451: EPSGetEigenvector - Gets the i-th right eigenvector as computed by EPSSolve().
453: Logically Collective on EPS455: Input Parameters:
456: + eps - eigensolver context
457: - i - index of the solution
459: Output Parameters:
460: + Vr - real part of eigenvector
461: - Vi - imaginary part of eigenvector
463: Notes:
464: If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
465: configured with complex scalars the eigenvector is stored
466: directly in Vr (Vi is set to zero).
468: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
469: Eigenpairs are indexed according to the ordering criterion established
470: with EPSSetWhichEigenpairs().
472: The 2-norm of the eigenvector is one unless the problem is generalized
473: Hermitian. In this case the eigenvector is normalized with respect to the
474: norm defined by the B matrix.
476: Level: beginner
478: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
479: @*/
480: PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi)481: {
483: PetscInt k;
491: EPSCheckSolved(eps,1);
492: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
493: EPSComputeVectors(eps);
494: if (!eps->perm) k = i;
495: else k = eps->perm[i];
496: #if defined(PETSC_USE_COMPLEX)
497: BVCopyVec(eps->V,k,Vr);
498: if (Vi) { VecSet(Vi,0.0); }
499: #else
500: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
501: BVCopyVec(eps->V,k,Vr);
502: if (Vi) {
503: BVCopyVec(eps->V,k+1,Vi);
504: }
505: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
506: BVCopyVec(eps->V,k-1,Vr);
507: if (Vi) {
508: BVCopyVec(eps->V,k,Vi);
509: VecScale(Vi,-1.0);
510: }
511: } else { /* real eigenvalue */
512: BVCopyVec(eps->V,k,Vr);
513: if (Vi) { VecSet(Vi,0.0); }
514: }
515: #endif
516: return(0);
517: }
521: /*@
522: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
523: computed eigenpair.
525: Not Collective
527: Input Parameter:
528: + eps - eigensolver context
529: - i - index of eigenpair
531: Output Parameter:
532: . errest - the error estimate
534: Notes:
535: This is the error estimate used internally by the eigensolver. The actual
536: error bound can be computed with EPSComputeRelativeError(). See also the users
537: manual for details.
539: Level: advanced
541: .seealso: EPSComputeRelativeError()
542: @*/
543: PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)544: {
548: EPSCheckSolved(eps,1);
549: if (i<0 || i>=eps->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
550: if (eps->perm) i = eps->perm[i];
551: if (errest) *errest = eps->errest[i];
552: return(0);
553: }
557: /*
558: EPSComputeResidualNorm_Private - Computes the norm of the residual vector
559: associated with an eigenpair.
560: */
561: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *norm)562: {
564: PetscInt nmat;
565: Vec u,w;
566: Mat A,B;
567: #if !defined(PETSC_USE_COMPLEX)
568: Vec v;
569: PetscReal ni,nr;
570: #endif
573: STGetNumMatrices(eps->st,&nmat);
574: STGetOperators(eps->st,0,&A);
575: if (nmat>1) { STGetOperators(eps->st,1,&B); }
576: BVGetVec(eps->V,&u);
577: BVGetVec(eps->V,&w);
579: #if !defined(PETSC_USE_COMPLEX)
580: if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
581: #endif
582: MatMult(A,xr,u); /* u=A*x */
583: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
584: if (eps->isgeneralized) { MatMult(B,xr,w); }
585: else { VecCopy(xr,w); } /* w=B*x */
586: VecAXPY(u,-kr,w); /* u=A*x-k*B*x */
587: }
588: VecNorm(u,NORM_2,norm);
589: #if !defined(PETSC_USE_COMPLEX)
590: } else {
591: BVGetVec(eps->V,&v);
592: MatMult(A,xr,u); /* u=A*xr */
593: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
594: if (eps->isgeneralized) { MatMult(B,xr,v); }
595: else { VecCopy(xr,v); } /* v=B*xr */
596: VecAXPY(u,-kr,v); /* u=A*xr-kr*B*xr */
597: if (eps->isgeneralized) { MatMult(B,xi,w); }
598: else { VecCopy(xi,w); } /* w=B*xi */
599: VecAXPY(u,ki,w); /* u=A*xr-kr*B*xr+ki*B*xi */
600: }
601: VecNorm(u,NORM_2,&nr);
602: MatMult(A,xi,u); /* u=A*xi */
603: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
604: VecAXPY(u,-kr,w); /* u=A*xi-kr*B*xi */
605: VecAXPY(u,-ki,v); /* u=A*xi-kr*B*xi-ki*B*xr */
606: }
607: VecNorm(u,NORM_2,&ni);
608: *norm = SlepcAbsEigenvalue(nr,ni);
609: VecDestroy(&v);
610: }
611: #endif
613: VecDestroy(&w);
614: VecDestroy(&u);
615: return(0);
616: }
620: /*@
621: EPSComputeResidualNorm - Computes the norm of the residual vector associated
622: with the i-th computed eigenpair.
624: Collective on EPS626: Input Parameter:
627: + eps - the eigensolver context
628: - i - the solution index
630: Output Parameter:
631: . norm - the residual norm, computed as ||Ax-kBx||_2 where k is the
632: eigenvalue and x is the eigenvector.
633: If k=0 then the residual norm is computed as ||Ax||_2.
635: Notes:
636: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
637: Eigenpairs are indexed according to the ordering criterion established
638: with EPSSetWhichEigenpairs().
640: Level: beginner
642: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
643: @*/
644: PetscErrorCode EPSComputeResidualNorm(EPS eps,PetscInt i,PetscReal *norm)645: {
647: Vec xr,xi;
648: PetscScalar kr,ki;
654: EPSCheckSolved(eps,1);
655: BVGetVec(eps->V,&xr);
656: BVGetVec(eps->V,&xi);
657: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
658: EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,norm);
659: VecDestroy(&xr);
660: VecDestroy(&xi);
661: return(0);
662: }
666: /*
667: EPSComputeRelativeError_Private - Computes the relative error bound
668: associated with an eigenpair.
669: */
670: PetscErrorCode EPSComputeRelativeError_Private(EPS eps,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,PetscReal *error)671: {
673: PetscReal norm,er;
674: #if !defined(PETSC_USE_COMPLEX)
675: PetscReal ei;
676: #endif
679: EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,&norm);
681: #if !defined(PETSC_USE_COMPLEX)
682: if (ki == 0) {
683: #endif
684: VecNorm(xr,NORM_2,&er);
685: #if !defined(PETSC_USE_COMPLEX)
686: } else {
687: VecNorm(xr,NORM_2,&er);
688: VecNorm(xi,NORM_2,&ei);
689: er = SlepcAbsEigenvalue(er,ei);
690: }
691: #endif
692: (*eps->converged)(eps,kr,ki,norm/er,error,eps->convergedctx);
693: return(0);
694: }
698: /*@
699: EPSComputeRelativeError - Computes the relative error bound associated
700: with the i-th computed eigenpair.
702: Collective on EPS704: Input Parameter:
705: + eps - the eigensolver context
706: - i - the solution index
708: Output Parameter:
709: . error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where
710: k is the eigenvalue and x is the eigenvector.
711: If k=0 the relative error is computed as ||Ax||_2/||x||_2.
713: Level: beginner
715: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
716: @*/
717: PetscErrorCode EPSComputeRelativeError(EPS eps,PetscInt i,PetscReal *error)718: {
720: Vec xr,xi;
721: PetscScalar kr,ki;
727: EPSCheckSolved(eps,1);
728: BVGetVec(eps->V,&xr);
729: BVGetVec(eps->V,&xi);
730: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
731: EPSComputeRelativeError_Private(eps,kr,ki,xr,xi,error);
732: VecDestroy(&xr);
733: VecDestroy(&xi);
734: return(0);
735: }
739: /*
740: EPSGetStartVector - Generate a suitable vector to be used as the starting vector
741: for the recurrence that builds the right subspace.
743: Collective on EPS and Vec
745: Input Parameters:
746: + eps - the eigensolver context
747: - i - iteration number
749: Output Parameters:
750: . breakdown - flag indicating that a breakdown has occurred
752: Notes:
753: The start vector is computed from another vector: for the first step (i=0),
754: the first initial vector is used (see EPSSetInitialSpace()); otherwise a random
755: vector is created. Then this vector is forced to be in the range of OP (only
756: for generalized definite problems) and orthonormalized with respect to all
757: V-vectors up to i-1. The resulting vector is placed in V[i].
759: The flag breakdown is set to true if either i=0 and the vector belongs to the
760: deflation space, or i>0 and the vector is linearly dependent with respect
761: to the V-vectors.
762: */
763: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,PetscBool *breakdown)764: {
766: PetscReal norm;
767: PetscBool lindep;
768: Vec w,z;
774: /* For the first step, use the first initial vector, otherwise a random one */
775: if (i>0 || eps->nini==0) {
776: BVSetRandomColumn(eps->V,i,eps->rand);
777: }
778: BVGetVec(eps->V,&w);
779: BVCopyVec(eps->V,i,w);
781: /* Force the vector to be in the range of OP for definite generalized problems */
782: BVGetColumn(eps->V,i,&z);
783: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
784: STApply(eps->st,w,z);
785: } else {
786: VecCopy(w,z);
787: }
788: BVRestoreColumn(eps->V,i,&z);
789: VecDestroy(&w);
791: /* Orthonormalize the vector with respect to previous vectors */
792: BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&lindep);
793: if (breakdown) *breakdown = lindep;
794: else if (lindep || norm == 0.0) {
795: if (i==0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Initial vector is zero or belongs to the deflation space");
796: else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to generate more start vectors");
797: }
798: BVScaleColumn(eps->V,i,1.0/norm);
799: return(0);
800: }