1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: EPS routines related to the solution process
12: */
14: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
15: #include <petscdraw.h>
17: PetscErrorCode EPSComputeVectors(EPS eps) 18: {
22: EPSCheckSolved(eps,1);
23: if (eps->state==EPS_STATE_SOLVED && eps->ops->computevectors) {
24: (*eps->ops->computevectors)(eps);
25: }
26: eps->state = EPS_STATE_EIGENVECTORS;
27: return(0);
28: }
30: #define SWAP(a,b,t) {t=a;a=b;b=t;} 32: static PetscErrorCode EPSComputeValues(EPS eps) 33: {
35: PetscBool injective,iscomp,isfilter;
36: PetscInt i,n,aux,nconv0;
37: Mat A,B=NULL,G,Z;
40: switch (eps->categ) {
41: case EPS_CATEGORY_KRYLOV:
42: case EPS_CATEGORY_OTHER:
43: STIsInjective(eps->st,&injective);
44: if (injective) {
45: /* one-to-one mapping: backtransform eigenvalues */
46: if (eps->ops->backtransform) {
47: (*eps->ops->backtransform)(eps);
48: } else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, spectral transform should have a backtransform operation");
49: } else {
50: /* compute eigenvalues from Rayleigh quotient */
51: DSGetDimensions(eps->ds,&n,NULL,NULL,NULL,NULL);
52: if (!n) break;
53: EPSGetOperators(eps,&A,&B);
54: BVSetActiveColumns(eps->V,0,n);
55: DSGetCompact(eps->ds,&iscomp);
56: DSSetCompact(eps->ds,PETSC_FALSE);
57: DSGetMat(eps->ds,DS_MAT_A,&G);
58: BVMatProject(eps->V,A,eps->V,G);
59: DSRestoreMat(eps->ds,DS_MAT_A,&G);
60: if (B) {
61: DSGetMat(eps->ds,DS_MAT_B,&G);
62: BVMatProject(eps->V,B,eps->V,G);
63: DSRestoreMat(eps->ds,DS_MAT_A,&G);
64: }
65: DSSolve(eps->ds,eps->eigr,eps->eigi);
66: DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
67: DSSynchronize(eps->ds,eps->eigr,eps->eigi);
68: DSSetCompact(eps->ds,iscomp);
69: if (eps->ishermitian && (!eps->isgeneralized || eps->ispositive)) { /* V = V * Z */
70: DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
71: DSGetMat(eps->ds,DS_MAT_X,&Z);
72: BVMultInPlace(eps->V,Z,0,n);
73: MatDestroy(&Z);
74: }
75: /* in case of STFILTER discard computed eigenvalues that lie outside the wanted interval */
76: PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilter);
77: if (isfilter) {
78: nconv0 = eps->nconv;
79: for (i=0;i<eps->nconv;i++) {
80: if (PetscRealPart(eps->eigr[eps->perm[i]])<eps->inta || PetscRealPart(eps->eigr[eps->perm[i]])>eps->intb) {
81: eps->nconv--;
82: if (i<eps->nconv) { SWAP(eps->perm[i],eps->perm[eps->nconv],aux); i--; }
83: }
84: }
85: if (nconv0>eps->nconv) {
86: PetscInfo1(eps,"Discarded %D computed eigenvalues lying outside the interval\n",nconv0-eps->nconv);
87: }
88: }
89: }
90: break;
91: case EPS_CATEGORY_PRECOND:
92: case EPS_CATEGORY_CONTOUR:
93: /* eigenvalues already available as an output of the solver */
94: break;
95: }
96: return(0);
97: }
99: /*@
100: EPSSolve - Solves the eigensystem.
102: Collective on EPS104: Input Parameter:
105: . eps - eigensolver context obtained from EPSCreate()
107: Options Database Keys:
108: + -eps_view - print information about the solver used
109: . -eps_view_mat0 binary - save the first matrix (A) to the default binary viewer
110: . -eps_view_mat1 binary - save the second matrix (B) to the default binary viewer
111: . -eps_view_vectors binary - save the computed eigenvectors to the default binary viewer
112: . -eps_view_values - print computed eigenvalues
113: . -eps_converged_reason - print reason for convergence, and number of iterations
114: . -eps_error_absolute - print absolute errors of each eigenpair
115: . -eps_error_relative - print relative errors of each eigenpair
116: - -eps_error_backward - print backward errors of each eigenpair
118: Level: beginner
120: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
121: @*/
122: PetscErrorCode EPSSolve(EPS eps)123: {
125: PetscInt i;
126: STMatMode matmode;
127: Mat A,B;
131: if (eps->state>=EPS_STATE_SOLVED) return(0);
132: PetscLogEventBegin(EPS_Solve,eps,0,0,0);
134: /* call setup */
135: EPSSetUp(eps);
136: eps->nconv = 0;
137: eps->its = 0;
138: for (i=0;i<eps->ncv;i++) {
139: eps->eigr[i] = 0.0;
140: eps->eigi[i] = 0.0;
141: eps->errest[i] = 0.0;
142: eps->perm[i] = i;
143: }
144: EPSViewFromOptions(eps,NULL,"-eps_view_pre");
146: /* call solver */
147: (*eps->ops->solve)(eps);
148: eps->state = EPS_STATE_SOLVED;
150: STGetMatMode(eps->st,&matmode);
151: if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
152: /* Purify eigenvectors before reverting operator */
153: EPSComputeVectors(eps);
154: }
155: STPostSolve(eps->st);
157: if (!eps->reason) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
159: /* Map eigenvalues back to the original problem if appropriate */
160: EPSComputeValues(eps);
162: #if !defined(PETSC_USE_COMPLEX)
163: /* reorder conjugate eigenvalues (positive imaginary first) */
164: for (i=0;i<eps->nconv-1;i++) {
165: if (eps->eigi[i] != 0) {
166: if (eps->eigi[i] < 0) {
167: eps->eigi[i] = -eps->eigi[i];
168: eps->eigi[i+1] = -eps->eigi[i+1];
169: /* the next correction only works with eigenvectors */
170: EPSComputeVectors(eps);
171: BVScaleColumn(eps->V,i+1,-1.0);
172: }
173: i++;
174: }
175: }
176: #endif
178: /* sort eigenvalues according to eps->which parameter */
179: SlepcSortEigenvalues(eps->sc,eps->nconv,eps->eigr,eps->eigi,eps->perm);
180: PetscLogEventEnd(EPS_Solve,eps,0,0,0);
182: /* various viewers */
183: EPSViewFromOptions(eps,NULL,"-eps_view");
184: EPSReasonViewFromOptions(eps);
185: EPSErrorViewFromOptions(eps);
186: EPSValuesViewFromOptions(eps);
187: EPSVectorsViewFromOptions(eps);
188: EPSGetOperators(eps,&A,&B);
189: MatViewFromOptions(A,(PetscObject)eps,"-eps_view_mat0");
190: if (eps->isgeneralized) {
191: MatViewFromOptions(B,(PetscObject)eps,"-eps_view_mat1");
192: }
194: /* Remove deflation and initial subspaces */
195: if (eps->nds) {
196: BVSetNumConstraints(eps->V,0);
197: eps->nds = 0;
198: }
199: eps->nini = 0;
200: return(0);
201: }
203: /*@
204: EPSGetIterationNumber - Gets the current iteration number. If the
205: call to EPSSolve() is complete, then it returns the number of iterations
206: carried out by the solution method.
208: Not Collective
210: Input Parameter:
211: . eps - the eigensolver context
213: Output Parameter:
214: . its - number of iterations
216: Note:
217: During the i-th iteration this call returns i-1. If EPSSolve() is
218: complete, then parameter "its" contains either the iteration number at
219: which convergence was successfully reached, or failure was detected.
220: Call EPSGetConvergedReason() to determine if the solver converged or
221: failed and why.
223: Level: intermediate
225: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
226: @*/
227: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)228: {
232: *its = eps->its;
233: return(0);
234: }
236: /*@
237: EPSGetConverged - Gets the number of converged eigenpairs.
239: Not Collective
241: Input Parameter:
242: . eps - the eigensolver context
244: Output Parameter:
245: . nconv - number of converged eigenpairs
247: Note:
248: This function should be called after EPSSolve() has finished.
250: Level: beginner
252: .seealso: EPSSetDimensions(), EPSSolve()
253: @*/
254: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)255: {
259: EPSCheckSolved(eps,1);
260: *nconv = eps->nconv;
261: return(0);
262: }
264: /*@
265: EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
266: stopped.
268: Not Collective
270: Input Parameter:
271: . eps - the eigensolver context
273: Output Parameter:
274: . reason - negative value indicates diverged, positive value converged
276: Notes:
277: Possible values for reason are
278: + EPS_CONVERGED_TOL - converged up to tolerance
279: . EPS_CONVERGED_USER - converged due to a user-defined condition
280: . EPS_DIVERGED_ITS - required more than max_it iterations to reach convergence
281: . EPS_DIVERGED_BREAKDOWN - generic breakdown in method
282: - EPS_DIVERGED_SYMMETRY_LOST - pseudo-Lanczos was not able to keep symmetry
284: Can only be called after the call to EPSSolve() is complete.
286: Level: intermediate
288: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason289: @*/
290: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)291: {
295: EPSCheckSolved(eps,1);
296: *reason = eps->reason;
297: return(0);
298: }
300: /*@
301: EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
302: subspace.
304: Not Collective, but vectors are shared by all processors that share the EPS306: Input Parameter:
307: . eps - the eigensolver context
309: Output Parameter:
310: . v - an array of vectors
312: Notes:
313: This function should be called after EPSSolve() has finished.
315: The user should provide in v an array of nconv vectors, where nconv is
316: the value returned by EPSGetConverged().
318: The first k vectors returned in v span an invariant subspace associated
319: with the first k computed eigenvalues (note that this is not true if the
320: k-th eigenvalue is complex and matrix A is real; in this case the first
321: k+1 vectors should be used). An invariant subspace X of A satisfies Ax
322: in X for all x in X (a similar definition applies for generalized
323: eigenproblems).
325: Level: intermediate
327: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve()
328: @*/
329: PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec *v)330: {
332: PetscInt i;
338: EPSCheckSolved(eps,1);
339: if (!eps->ishermitian && eps->state==EPS_STATE_EIGENVECTORS) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector or EPSComputeError");
340: for (i=0;i<eps->nconv;i++) {
341: BVCopyVec(eps->V,i,v[i]);
342: if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
343: VecPointwiseDivide(v[i],v[i],eps->D);
344: VecNormalize(v[i],NULL);
345: }
346: }
347: return(0);
348: }
350: /*@C
351: EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
352: EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.
354: Logically Collective on EPS356: Input Parameters:
357: + eps - eigensolver context
358: - i - index of the solution
360: Output Parameters:
361: + eigr - real part of eigenvalue
362: . eigi - imaginary part of eigenvalue
363: . Vr - real part of eigenvector
364: - Vi - imaginary part of eigenvector
366: Notes:
367: It is allowed to pass NULL for Vr and Vi, if the eigenvector is not
368: required. Otherwise, the caller must provide valid Vec objects, i.e.,
369: they must be created by the calling program with e.g. MatCreateVecs().
371: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
372: configured with complex scalars the eigenvalue is stored
373: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
374: set to zero). In both cases, the user can pass NULL in eigi and Vi.
376: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
377: Eigenpairs are indexed according to the ordering criterion established
378: with EPSSetWhichEigenpairs().
380: The 2-norm of the eigenvector is one unless the problem is generalized
381: Hermitian. In this case the eigenvector is normalized with respect to the
382: norm defined by the B matrix.
384: Level: beginner
386: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSSolve(),
387: EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
388: @*/
389: PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)390: {
396: EPSCheckSolved(eps,1);
397: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
398: EPSGetEigenvalue(eps,i,eigr,eigi);
399: if (Vr || Vi) { EPSGetEigenvector(eps,i,Vr,Vi); }
400: return(0);
401: }
403: /*@C
404: EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve().
406: Not Collective
408: Input Parameters:
409: + eps - eigensolver context
410: - i - index of the solution
412: Output Parameters:
413: + eigr - real part of eigenvalue
414: - eigi - imaginary part of eigenvalue
416: Notes:
417: If the eigenvalue is real, then eigi is set to zero. If PETSc is
418: configured with complex scalars the eigenvalue is stored
419: directly in eigr (eigi is set to zero).
421: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
422: Eigenpairs are indexed according to the ordering criterion established
423: with EPSSetWhichEigenpairs().
425: Level: beginner
427: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
428: @*/
429: PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi)430: {
431: PetscInt k;
435: EPSCheckSolved(eps,1);
436: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
437: 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: }
448: /*@C
449: EPSGetEigenvector - Gets the i-th right eigenvector as computed by EPSSolve().
451: Logically Collective on EPS453: Input Parameters:
454: + eps - eigensolver context
455: - i - index of the solution
457: Output Parameters:
458: + Vr - real part of eigenvector
459: - Vi - imaginary part of eigenvector
461: Notes:
462: The caller must provide valid Vec objects, i.e., they must be created
463: by the calling program with e.g. MatCreateVecs().
465: If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
466: configured with complex scalars the eigenvector is stored
467: directly in Vr (Vi is set to zero). In both cases, the user can pass NULL in Vi.
469: The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
470: Eigenpairs are indexed according to the ordering criterion established
471: with EPSSetWhichEigenpairs().
473: The 2-norm of the eigenvector is one unless the problem is generalized
474: Hermitian. In this case the eigenvector is normalized with respect to the
475: norm defined by the B matrix.
477: Level: beginner
479: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
480: @*/
481: PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi)482: {
484: PetscInt k;
492: EPSCheckSolved(eps,1);
493: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
494: EPSComputeVectors(eps);
495: 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: }
519: /*@
520: EPSGetErrorEstimate - Returns the error estimate associated to the i-th
521: computed eigenpair.
523: Not Collective
525: Input Parameter:
526: + eps - eigensolver context
527: - i - index of eigenpair
529: Output Parameter:
530: . errest - the error estimate
532: Notes:
533: This is the error estimate used internally by the eigensolver. The actual
534: error bound can be computed with EPSComputeError(). See also the users
535: manual for details.
537: Level: advanced
539: .seealso: EPSComputeError()
540: @*/
541: PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)542: {
546: EPSCheckSolved(eps,1);
547: if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
548: *errest = eps->errest[eps->perm[i]];
549: return(0);
550: }
552: /*
553: EPSComputeResidualNorm_Private - Computes the norm of the residual vector
554: associated with an eigenpair.
556: Input Parameters:
557: kr,ki - eigenvalue
558: xr,xi - eigenvector
559: z - three work vectors (the second one not referenced in complex scalars)
560: */
561: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)562: {
564: PetscInt nmat;
565: Mat A,B;
566: Vec u,w;
567: #if !defined(PETSC_USE_COMPLEX)
568: Vec v;
569: PetscReal ni,nr;
570: #endif
573: u = z[0]; w = z[2];
574: STGetNumMatrices(eps->st,&nmat);
575: STGetMatrix(eps->st,0,&A);
576: if (nmat>1) { STGetMatrix(eps->st,1,&B); }
578: #if !defined(PETSC_USE_COMPLEX)
579: v = z[1];
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 (nmat>1) { 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: MatMult(A,xr,u); /* u=A*xr */
592: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
593: if (nmat>1) { MatMult(B,xr,v); }
594: else { VecCopy(xr,v); } /* v=B*xr */
595: VecAXPY(u,-kr,v); /* u=A*xr-kr*B*xr */
596: if (nmat>1) { MatMult(B,xi,w); }
597: else { VecCopy(xi,w); } /* w=B*xi */
598: VecAXPY(u,ki,w); /* u=A*xr-kr*B*xr+ki*B*xi */
599: }
600: VecNorm(u,NORM_2,&nr);
601: MatMult(A,xi,u); /* u=A*xi */
602: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
603: VecAXPY(u,-kr,w); /* u=A*xi-kr*B*xi */
604: VecAXPY(u,-ki,v); /* u=A*xi-kr*B*xi-ki*B*xr */
605: }
606: VecNorm(u,NORM_2,&ni);
607: *norm = SlepcAbsEigenvalue(nr,ni);
608: }
609: #endif
610: return(0);
611: }
613: /*@
614: EPSComputeError - Computes the error (based on the residual norm) associated
615: with the i-th computed eigenpair.
617: Collective on EPS619: Input Parameter:
620: + eps - the eigensolver context
621: . i - the solution index
622: - type - the type of error to compute
624: Output Parameter:
625: . error - the error
627: Notes:
628: The error can be computed in various ways, all of them based on the residual
629: norm ||Ax-kBx||_2 where k is the eigenvalue and x is the eigenvector.
631: Level: beginner
633: .seealso: EPSErrorType, EPSSolve(), EPSGetErrorEstimate()
634: @*/
635: PetscErrorCode EPSComputeError(EPS eps,PetscInt i,EPSErrorType type,PetscReal *error)636: {
638: Mat A,B;
639: Vec xr,xi,w[3];
640: PetscReal t,vecnorm=1.0;
641: PetscScalar kr,ki;
642: PetscBool flg;
649: EPSCheckSolved(eps,1);
651: /* allocate work vectors */
652: #if defined(PETSC_USE_COMPLEX)
653: EPSSetWorkVecs(eps,3);
654: xi = NULL;
655: w[1] = NULL;
656: #else
657: EPSSetWorkVecs(eps,5);
658: xi = eps->work[3];
659: w[1] = eps->work[4];
660: #endif
661: xr = eps->work[0];
662: w[0] = eps->work[1];
663: w[2] = eps->work[2];
665: /* compute residual norms */
666: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
667: EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,w,error);
669: /* compute 2-norm of eigenvector */
670: if (eps->problem_type==EPS_GHEP) {
671: VecNorm(xr,NORM_2,&vecnorm);
672: }
674: /* compute error */
675: switch (type) {
676: case EPS_ERROR_ABSOLUTE:
677: break;
678: case EPS_ERROR_RELATIVE:
679: *error /= SlepcAbsEigenvalue(kr,ki)*vecnorm;
680: break;
681: case EPS_ERROR_BACKWARD:
682: /* initialization of matrix norms */
683: if (!eps->nrma) {
684: STGetMatrix(eps->st,0,&A);
685: MatHasOperation(A,MATOP_NORM,&flg);
686: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
687: MatNorm(A,NORM_INFINITY,&eps->nrma);
688: }
689: if (eps->isgeneralized) {
690: if (!eps->nrmb) {
691: STGetMatrix(eps->st,1,&B);
692: MatHasOperation(B,MATOP_NORM,&flg);
693: if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
694: MatNorm(B,NORM_INFINITY,&eps->nrmb);
695: }
696: } else eps->nrmb = 1.0;
697: t = SlepcAbsEigenvalue(kr,ki);
698: *error /= (eps->nrma+t*eps->nrmb)*vecnorm;
699: break;
700: default:701: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
702: }
703: return(0);
704: }
706: /*
707: EPSGetStartVector - Generate a suitable vector to be used as the starting vector
708: for the recurrence that builds the right subspace.
710: Collective on EPS and Vec
712: Input Parameters:
713: + eps - the eigensolver context
714: - i - iteration number
716: Output Parameters:
717: . breakdown - flag indicating that a breakdown has occurred
719: Notes:
720: The start vector is computed from another vector: for the first step (i=0),
721: the first initial vector is used (see EPSSetInitialSpace()); otherwise a random
722: vector is created. Then this vector is forced to be in the range of OP (only
723: for generalized definite problems) and orthonormalized with respect to all
724: V-vectors up to i-1. The resulting vector is placed in V[i].
726: The flag breakdown is set to true if either i=0 and the vector belongs to the
727: deflation space, or i>0 and the vector is linearly dependent with respect
728: to the V-vectors.
729: */
730: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,PetscBool *breakdown)731: {
733: PetscReal norm;
734: PetscBool lindep;
735: Vec w,z;
741: /* For the first step, use the first initial vector, otherwise a random one */
742: if (i>0 || eps->nini==0) {
743: BVSetRandomColumn(eps->V,i);
744: }
746: /* Force the vector to be in the range of OP for definite generalized problems */
747: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
748: BVCreateVec(eps->V,&w);
749: BVCopyVec(eps->V,i,w);
750: BVGetColumn(eps->V,i,&z);
751: STApply(eps->st,w,z);
752: BVRestoreColumn(eps->V,i,&z);
753: VecDestroy(&w);
754: }
756: /* Orthonormalize the vector with respect to previous vectors */
757: BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&lindep);
758: if (breakdown) *breakdown = lindep;
759: else if (lindep || norm == 0.0) {
760: if (i==0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Initial vector is zero or belongs to the deflation space");
761: else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to generate more start vectors");
762: }
763: BVScaleColumn(eps->V,i,1.0/norm);
764: return(0);
765: }