Actual source code: epssolve.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:       EPS routines related to the solution process.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, 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: 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_view_vectors binary - save the computed eigenvectors to the default binary viewer
 63: .  -eps_view_values - print computed eigenvalues
 64: .  -eps_converged_reason - print reason for convergence, and number of iterations
 65: .  -eps_error_absolute - print absolute errors of each eigenpair
 66: .  -eps_error_relative - print relative errors of each eigenpair
 67: -  -eps_error_backward - print backward errors of each eigenpair

 69:    Level: beginner

 71: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
 72: @*/
 73: PetscErrorCode EPSSolve(EPS eps)
 74: {
 76:   PetscInt       i,nmat;
 77:   PetscScalar    dot;
 78:   PetscBool      iscayley;
 79:   STMatMode      matmode;
 80:   Mat            A,B;
 81:   Vec            w,x;

 85:   if (eps->state>=EPS_STATE_SOLVED) return(0);
 86:   PetscLogEventBegin(EPS_Solve,eps,0,0,0);

 88:   /* call setup */
 89:   EPSSetUp(eps);
 90:   eps->nconv = 0;
 91:   eps->its   = 0;
 92:   for (i=0;i<eps->ncv;i++) {
 93:     eps->eigr[i]   = 0.0;
 94:     eps->eigi[i]   = 0.0;
 95:     eps->errest[i] = 0.0;
 96:     eps->perm[i]   = i;
 97:   }
 98:   EPSViewFromOptions(eps,NULL,"-eps_view_pre");

100:   /* call solver */
101:   (*eps->ops->solve)(eps);
102:   eps->state = EPS_STATE_SOLVED;

104:   STGetMatMode(eps->st,&matmode);
105:   if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
106:     /* Purify eigenvectors before reverting operator */
107:     EPSComputeVectors(eps);
108:   }
109:   STPostSolve(eps->st);

111:   if (!eps->reason) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");

113:   /* Map eigenvalues back to the original problem, necessary in some
114:   * spectral transformations */
115:   if (eps->ops->backtransform) {
116:     (*eps->ops->backtransform)(eps);
117:   }

119: #if !defined(PETSC_USE_COMPLEX)
120:   /* reorder conjugate eigenvalues (positive imaginary first) */
121:   for (i=0; i<eps->nconv-1; i++) {
122:     if (eps->eigi[i] != 0) {
123:       if (eps->eigi[i] < 0) {
124:         eps->eigi[i] = -eps->eigi[i];
125:         eps->eigi[i+1] = -eps->eigi[i+1];
126:         /* the next correction only works with eigenvectors */
127:         EPSComputeVectors(eps);
128:         BVScaleColumn(eps->V,i+1,-1.0);
129:       }
130:       i++;
131:     }
132:   }
133: #endif

135:   STGetNumMatrices(eps->st,&nmat);
136:   STGetOperators(eps->st,0,&A);
137:   if (nmat>1) { STGetOperators(eps->st,1,&B); }

139:   /* In the case of Cayley transform, eigenvectors need to be B-normalized */
140:   PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&iscayley);
141:   if (iscayley && nmat>1 && eps->ishermitian) {
142:     MatCreateVecs(B,NULL,&w);
143:     EPSComputeVectors(eps);
144:     for (i=0;i<eps->nconv;i++) {
145:       BVGetColumn(eps->V,i,&x);
146:       MatMult(B,x,w);
147:       VecDot(w,x,&dot);
148:       VecScale(x,1.0/PetscSqrtScalar(dot));
149:       BVRestoreColumn(eps->V,i,&x);
150:     }
151:     VecDestroy(&w);
152:   }

154:   /* sort eigenvalues according to eps->which parameter */
155:   SlepcSortEigenvalues(eps->sc,eps->nconv,eps->eigr,eps->eigi,eps->perm);
156:   PetscLogEventEnd(EPS_Solve,eps,0,0,0);

158:   /* various viewers */
159:   EPSViewFromOptions(eps,NULL,"-eps_view");
160:   EPSReasonViewFromOptions(eps);
161:   EPSErrorViewFromOptions(eps);
162:   EPSValuesViewFromOptions(eps);
163:   EPSVectorsViewFromOptions(eps);
164:   MatViewFromOptions(A,(PetscObject)eps,"-eps_view_mat0");
165:   if (nmat>1) { MatViewFromOptions(B,(PetscObject)eps,"-eps_view_mat1"); }

167:   /* Remove deflation and initial subspaces */
168:   if (eps->nds) {
169:     BVSetNumConstraints(eps->V,0);
170:     eps->nds = 0;
171:   }
172:   eps->nini = 0;
173:   return(0);
174: }

178: /*@
179:    EPSGetIterationNumber - Gets the current iteration number. If the
180:    call to EPSSolve() is complete, then it returns the number of iterations
181:    carried out by the solution method.

183:    Not Collective

185:    Input Parameter:
186: .  eps - the eigensolver context

188:    Output Parameter:
189: .  its - number of iterations

191:    Level: intermediate

193:    Note:
194:    During the i-th iteration this call returns i-1. If EPSSolve() is
195:    complete, then parameter "its" contains either the iteration number at
196:    which convergence was successfully reached, or failure was detected.
197:    Call EPSGetConvergedReason() to determine if the solver converged or
198:    failed and why.

200: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
201: @*/
202: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
203: {
207:   *its = eps->its;
208:   return(0);
209: }

213: /*@
214:    EPSGetConverged - Gets the number of converged eigenpairs.

216:    Not Collective

218:    Input Parameter:
219: .  eps - the eigensolver context

221:    Output Parameter:
222: .  nconv - number of converged eigenpairs

224:    Note:
225:    This function should be called after EPSSolve() has finished.

227:    Level: beginner

229: .seealso: EPSSetDimensions(), EPSSolve()
230: @*/
231: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
232: {
236:   EPSCheckSolved(eps,1);
237:   *nconv = eps->nconv;
238:   return(0);
239: }

243: /*@
244:    EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
245:    stopped.

247:    Not Collective

249:    Input Parameter:
250: .  eps - the eigensolver context

252:    Output Parameter:
253: .  reason - negative value indicates diverged, positive value converged

255:    Possible values for reason:
256: +  EPS_CONVERGED_TOL - converged up to tolerance
257: .  EPS_CONVERGED_USER - converged due to a user-defined condition
258: .  EPS_DIVERGED_ITS - required more than max_it iterations to reach convergence
259: .  EPS_DIVERGED_BREAKDOWN - generic breakdown in method
260: -  EPS_DIVERGED_SYMMETRY_LOST - pseudo-Lanczos was not able to keep symmetry

262:    Note:
263:    Can only be called after the call to EPSSolve() is complete.

265:    Level: intermediate

267: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
268: @*/
269: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
270: {
274:   EPSCheckSolved(eps,1);
275:   *reason = eps->reason;
276:   return(0);
277: }

281: /*@
282:    EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
283:    subspace.

285:    Not Collective, but vectors are shared by all processors that share the EPS

287:    Input Parameter:
288: .  eps - the eigensolver context

290:    Output Parameter:
291: .  v - an array of vectors

293:    Notes:
294:    This function should be called after EPSSolve() has finished.

296:    The user should provide in v an array of nconv vectors, where nconv is
297:    the value returned by EPSGetConverged().

299:    The first k vectors returned in v span an invariant subspace associated
300:    with the first k computed eigenvalues (note that this is not true if the
301:    k-th eigenvalue is complex and matrix A is real; in this case the first
302:    k+1 vectors should be used). An invariant subspace X of A satisfies Ax
303:    in X for all x in X (a similar definition applies for generalized
304:    eigenproblems).

306:    Level: intermediate

308: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve()
309: @*/
310: PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec *v)
311: {
313:   PetscInt       i;

319:   EPSCheckSolved(eps,1);
320:   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");
321:   for (i=0;i<eps->nconv;i++) {
322:     BVCopyVec(eps->V,i,v[i]);
323:     if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
324:       VecPointwiseDivide(v[i],v[i],eps->D);
325:       VecNormalize(v[i],NULL);
326:     }
327:   }
328:   return(0);
329: }

333: /*@
334:    EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
335:    EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.

337:    Logically Collective on EPS

339:    Input Parameters:
340: +  eps - eigensolver context
341: -  i   - index of the solution

343:    Output Parameters:
344: +  eigr - real part of eigenvalue
345: .  eigi - imaginary part of eigenvalue
346: .  Vr   - real part of eigenvector
347: -  Vi   - imaginary part of eigenvector

349:    Notes:
350:    It is allowed to pass NULL for Vr and Vi, if the eigenvector is not
351:    required. Otherwise, the caller must provide valid Vec objects, i.e.,
352:    they must be created by the calling program with e.g. MatCreateVecs().

354:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
355:    configured with complex scalars the eigenvalue is stored
356:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
357:    set to zero). In both cases, the user can pass NULL in eigi and Vi.

359:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
360:    Eigenpairs are indexed according to the ordering criterion established
361:    with EPSSetWhichEigenpairs().

363:    The 2-norm of the eigenvector is one unless the problem is generalized
364:    Hermitian. In this case the eigenvector is normalized with respect to the
365:    norm defined by the B matrix.

367:    Level: beginner

369: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSSolve(),
370:           EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
371: @*/
372: PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
373: {

379:   EPSCheckSolved(eps,1);
380:   if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
381:   EPSGetEigenvalue(eps,i,eigr,eigi);
382:   if (Vr || Vi) { EPSGetEigenvector(eps,i,Vr,Vi); }
383:   return(0);
384: }

388: /*@
389:    EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve().

391:    Not Collective

393:    Input Parameters:
394: +  eps - eigensolver context
395: -  i   - index of the solution

397:    Output Parameters:
398: +  eigr - real part of eigenvalue
399: -  eigi - imaginary part of eigenvalue

401:    Notes:
402:    If the eigenvalue is real, then eigi is set to zero. If PETSc is
403:    configured with complex scalars the eigenvalue is stored
404:    directly in eigr (eigi is set to zero).

406:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
407:    Eigenpairs are indexed according to the ordering criterion established
408:    with EPSSetWhichEigenpairs().

410:    Level: beginner

412: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
413: @*/
414: PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi)
415: {
416:   PetscInt k;

420:   EPSCheckSolved(eps,1);
421:   if (i<0 || i>=eps->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
422:   k = eps->perm[i];
423: #if defined(PETSC_USE_COMPLEX)
424:   if (eigr) *eigr = eps->eigr[k];
425:   if (eigi) *eigi = 0;
426: #else
427:   if (eigr) *eigr = eps->eigr[k];
428:   if (eigi) *eigi = eps->eigi[k];
429: #endif
430:   return(0);
431: }

435: /*@
436:    EPSGetEigenvector - Gets the i-th right eigenvector as computed by EPSSolve().

438:    Logically Collective on EPS

440:    Input Parameters:
441: +  eps - eigensolver context
442: -  i   - index of the solution

444:    Output Parameters:
445: +  Vr   - real part of eigenvector
446: -  Vi   - imaginary part of eigenvector

448:    Notes:
449:    The caller must provide valid Vec objects, i.e., they must be created
450:    by the calling program with e.g. MatCreateVecs().

452:    If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
453:    configured with complex scalars the eigenvector is stored
454:    directly in Vr (Vi is set to zero). In both cases, the user can pass NULL in Vi.

456:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
457:    Eigenpairs are indexed according to the ordering criterion established
458:    with EPSSetWhichEigenpairs().

460:    The 2-norm of the eigenvector is one unless the problem is generalized
461:    Hermitian. In this case the eigenvector is normalized with respect to the
462:    norm defined by the B matrix.

464:    Level: beginner

466: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
467: @*/
468: PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi)
469: {
471:   PetscInt       k;

479:   EPSCheckSolved(eps,1);
480:   if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
481:   EPSComputeVectors(eps);
482:   k = eps->perm[i];
483: #if defined(PETSC_USE_COMPLEX)
484:   BVCopyVec(eps->V,k,Vr);
485:   if (Vi) { VecSet(Vi,0.0); }
486: #else
487:   if (eps->eigi[k] > 0) { /* first value of conjugate pair */
488:     BVCopyVec(eps->V,k,Vr);
489:     if (Vi) {
490:       BVCopyVec(eps->V,k+1,Vi);
491:     }
492:   } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
493:     BVCopyVec(eps->V,k-1,Vr);
494:     if (Vi) {
495:       BVCopyVec(eps->V,k,Vi);
496:       VecScale(Vi,-1.0);
497:     }
498:   } else { /* real eigenvalue */
499:     BVCopyVec(eps->V,k,Vr);
500:     if (Vi) { VecSet(Vi,0.0); }
501:   }
502: #endif
503:   return(0);
504: }

508: /*@
509:    EPSGetErrorEstimate - Returns the error estimate associated to the i-th
510:    computed eigenpair.

512:    Not Collective

514:    Input Parameter:
515: +  eps - eigensolver context
516: -  i   - index of eigenpair

518:    Output Parameter:
519: .  errest - the error estimate

521:    Notes:
522:    This is the error estimate used internally by the eigensolver. The actual
523:    error bound can be computed with EPSComputeError(). See also the users
524:    manual for details.

526:    Level: advanced

528: .seealso: EPSComputeError()
529: @*/
530: PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)
531: {
535:   EPSCheckSolved(eps,1);
536:   if (i<0 || i>=eps->nconv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
537:   if (errest) *errest = eps->errest[eps->perm[i]];
538:   return(0);
539: }

543: /*
544:    EPSComputeResidualNorm_Private - Computes the norm of the residual vector
545:    associated with an eigenpair.

547:    Input Parameters:
548:      kr,ki - eigenvalue
549:      xr,xi - eigenvector
550:      z     - three work vectors (the second one not referenced in complex scalars)
551: */
552: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
553: {
555:   PetscInt       nmat;
556:   Mat            A,B;
557:   Vec            u,w;
558: #if !defined(PETSC_USE_COMPLEX)
559:   Vec            v;
560:   PetscReal      ni,nr;
561: #endif

564:   u = z[0]; w = z[2];
565:   STGetNumMatrices(eps->st,&nmat);
566:   STGetOperators(eps->st,0,&A);
567:   if (nmat>1) { STGetOperators(eps->st,1,&B); }

569: #if !defined(PETSC_USE_COMPLEX)
570:   v = z[1]; 
571:   if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
572: #endif
573:     MatMult(A,xr,u);                             /* u=A*x */
574:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
575:       if (nmat>1) { MatMult(B,xr,w); }
576:       else { VecCopy(xr,w); }                    /* w=B*x */
577:       VecAXPY(u,-kr,w);                          /* u=A*x-k*B*x */
578:     }
579:     VecNorm(u,NORM_2,norm);
580: #if !defined(PETSC_USE_COMPLEX)
581:   } else {
582:     MatMult(A,xr,u);                             /* u=A*xr */
583:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
584:       if (nmat>1) { MatMult(B,xr,v); }
585:       else { VecCopy(xr,v); }                    /* v=B*xr */
586:       VecAXPY(u,-kr,v);                          /* u=A*xr-kr*B*xr */
587:       if (nmat>1) { MatMult(B,xi,w); }
588:       else { VecCopy(xi,w); }                    /* w=B*xi */
589:       VecAXPY(u,ki,w);                           /* u=A*xr-kr*B*xr+ki*B*xi */
590:     }
591:     VecNorm(u,NORM_2,&nr);
592:     MatMult(A,xi,u);                             /* u=A*xi */
593:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
594:       VecAXPY(u,-kr,w);                          /* u=A*xi-kr*B*xi */
595:       VecAXPY(u,-ki,v);                          /* u=A*xi-kr*B*xi-ki*B*xr */
596:     }
597:     VecNorm(u,NORM_2,&ni);
598:     *norm = SlepcAbsEigenvalue(nr,ni);
599:   }
600: #endif
601:   return(0);
602: }

606: /*@
607:    EPSComputeError - Computes the error (based on the residual norm) associated
608:    with the i-th computed eigenpair.

610:    Collective on EPS

612:    Input Parameter:
613: +  eps  - the eigensolver context
614: .  i    - the solution index
615: -  type - the type of error to compute

617:    Output Parameter:
618: .  error - the error

620:    Notes:
621:    The error can be computed in various ways, all of them based on the residual
622:    norm ||Ax-kBx||_2 where k is the eigenvalue and x is the eigenvector.

624:    Level: beginner

626: .seealso: EPSErrorType, EPSSolve(), EPSGetErrorEstimate()
627: @*/
628: PetscErrorCode EPSComputeError(EPS eps,PetscInt i,EPSErrorType type,PetscReal *error)
629: {
631:   Mat            A,B;
632:   Vec            xr,xi,w[3];
633:   PetscReal      t;
634:   PetscScalar    kr,ki;
635:   PetscBool      flg;

642:   EPSCheckSolved(eps,1);

644:   /* allocate work vectors */
645: #if defined(PETSC_USE_COMPLEX)
646:   EPSSetWorkVecs(eps,3);
647:   xi   = NULL;
648:   w[1] = NULL;
649: #else
650:   EPSSetWorkVecs(eps,5);
651:   xi   = eps->work[3];
652:   w[1] = eps->work[4];
653: #endif
654:   xr   = eps->work[0];
655:   w[0] = eps->work[1];
656:   w[2] = eps->work[2];

658:   /* compute residual norms */
659:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
660:   EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,w,error);

662:   /* compute error */
663:   switch (type) {
664:     case EPS_ERROR_ABSOLUTE:
665:       break;
666:     case EPS_ERROR_RELATIVE:
667:       *error /= SlepcAbsEigenvalue(kr,ki);
668:       break;
669:     case EPS_ERROR_BACKWARD:
670:       /* initialization of matrix norms */
671:       if (!eps->nrma) {
672:         STGetOperators(eps->st,0,&A);
673:         MatHasOperation(A,MATOP_NORM,&flg);
674:         if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
675:         MatNorm(A,NORM_INFINITY,&eps->nrma);
676:       }
677:       if (eps->isgeneralized) {
678:         if (!eps->nrmb) {
679:           STGetOperators(eps->st,1,&B);
680:           MatHasOperation(B,MATOP_NORM,&flg);
681:           if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
682:           MatNorm(B,NORM_INFINITY,&eps->nrmb);
683:         }
684:       } else eps->nrmb = 1.0;
685:       t = SlepcAbsEigenvalue(kr,ki);
686:       *error /= eps->nrma+t*eps->nrmb;
687:       break;
688:     default:
689:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
690:   }
691:   return(0);
692: }

696: /*
697:    EPSGetStartVector - Generate a suitable vector to be used as the starting vector
698:    for the recurrence that builds the right subspace.

700:    Collective on EPS and Vec

702:    Input Parameters:
703: +  eps - the eigensolver context
704: -  i   - iteration number

706:    Output Parameters:
707: .  breakdown - flag indicating that a breakdown has occurred

709:    Notes:
710:    The start vector is computed from another vector: for the first step (i=0),
711:    the first initial vector is used (see EPSSetInitialSpace()); otherwise a random
712:    vector is created. Then this vector is forced to be in the range of OP (only
713:    for generalized definite problems) and orthonormalized with respect to all
714:    V-vectors up to i-1. The resulting vector is placed in V[i].

716:    The flag breakdown is set to true if either i=0 and the vector belongs to the
717:    deflation space, or i>0 and the vector is linearly dependent with respect
718:    to the V-vectors.
719: */
720: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,PetscBool *breakdown)
721: {
723:   PetscReal      norm;
724:   PetscBool      lindep;
725:   Vec            w,z;


731:   /* For the first step, use the first initial vector, otherwise a random one */
732:   if (i>0 || eps->nini==0) {
733:     BVSetRandomColumn(eps->V,i);
734:   }
735:   BVCreateVec(eps->V,&w);
736:   BVCopyVec(eps->V,i,w);

738:   /* Force the vector to be in the range of OP for definite generalized problems */
739:   BVGetColumn(eps->V,i,&z);
740:   if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
741:     STApply(eps->st,w,z);
742:   } else {
743:     VecCopy(w,z);
744:   }
745:   BVRestoreColumn(eps->V,i,&z);
746:   VecDestroy(&w);

748:   /* Orthonormalize the vector with respect to previous vectors */
749:   BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&lindep);
750:   if (breakdown) *breakdown = lindep;
751:   else if (lindep || norm == 0.0) {
752:     if (i==0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Initial vector is zero or belongs to the deflation space");
753:     else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to generate more start vectors");
754:   }
755:   BVScaleColumn(eps->V,i,1.0/norm);
756:   return(0);
757: }