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-2009, Universidad Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24:  #include private/epsimpl.h

 28: /*@
 29:    EPSSolve - Solves the eigensystem.

 31:    Collective on EPS

 33:    Input Parameter:
 34: .  eps - eigensolver context obtained from EPSCreate()

 36:    Options Database:
 37: +   -eps_view - print information about the solver used
 38: .   -eps_view_binary - save the matrices to the default binary file
 39: -   -eps_plot_eigs - plot computed eigenvalues

 41:    Level: beginner

 43: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances() 
 44: @*/
 45: PetscErrorCode EPSSolve(EPS eps)
 46: {
 48:   PetscInt       i;
 49:   PetscReal      re,im;
 50:   PetscTruth     flg;
 51:   PetscViewer    viewer;
 52:   PetscDraw      draw;
 53:   PetscDrawSP    drawsp;
 54:   STMatMode      matmode;


 59:   PetscOptionsHasName(((PetscObject)eps)->prefix,"-eps_view_binary",&flg);
 60:   if (flg) {
 61:     Mat A,B;
 62:     STGetOperators(eps->OP,&A,&B);
 63:     MatView(A,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
 64:     if (B) MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
 65:   }

 67:   /* reset the convergence flag from the previous solves */
 68:   eps->reason = EPS_CONVERGED_ITERATING;

 70:   if (!eps->setupcalled){ EPSSetUp(eps); }
 71:   STResetOperationCounters(eps->OP);
 72:   IPResetOperationCounters(eps->ip);
 73:   eps->evecsavailable = PETSC_FALSE;
 74:   eps->nconv = 0;
 75:   eps->its = 0;
 76:   for (i=0;i<eps->ncv;i++) eps->eigr[i]=eps->eigi[i]=eps->errest[i]=0.0;
 77:   EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->ncv);

 79:   PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);

 81:   switch (eps->solverclass) {
 82:     case EPS_ONE_SIDE:
 83:       (*eps->ops->solve)(eps); break;
 84:     case EPS_TWO_SIDE:
 85:       if (eps->ops->solvets) {
 86:         (*eps->ops->solvets)(eps); break;
 87:       } else {
 88:         SETERRQ(1,"Two-sided version unavailable for this solver");
 89:       }
 90:     default:
 91:       SETERRQ(1,"Wrong value of eps->solverclass");
 92:   }

 94:   STGetMatMode(eps->OP,&matmode);
 95:   if (matmode == STMATMODE_INPLACE && eps->ispositive) {
 96:     /* Purge eigenvectors before reverting operator */
 97:     (*eps->ops->computevectors)(eps);
 98:   }
 99:   STPostSolve(eps->OP);
100:   PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);

102:   if (!eps->reason) {
103:     SETERRQ(1,"Internal error, solver returned without setting converged reason");
104:   }

106:   /* Map eigenvalues back to the original problem, necessary in some 
107:   * spectral transformations */
108:   (*eps->ops->backtransform)(eps);

110:   /* Adjust left eigenvectors in generalized problems: y = B^T y */
111:   if (eps->isgeneralized && eps->solverclass == EPS_TWO_SIDE) {
112:     Mat B;
113:     KSP ksp;
114:     Vec w;
115:     STGetOperators(eps->OP,PETSC_NULL,&B);
116:     KSPCreate(((PetscObject)eps)->comm,&ksp);
117:     KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
118:     KSPSetFromOptions(ksp);
119:     MatGetVecs(B,PETSC_NULL,&w);
120:     for (i=0;i<eps->nconv;i++) {
121:       VecCopy(eps->W[i],w);
122:       KSPSolveTranspose(ksp,w,eps->W[i]);
123:     }
124:     KSPDestroy(ksp);
125:     VecDestroy(w);
126:   }
127: 

129: #ifndef PETSC_USE_COMPLEX
130:   /* reorder conjugate eigenvalues (positive imaginary first) */
131:   for (i=0; i<eps->nconv-1; i++) {
132:     if (eps->eigi[i] != 0) {
133:       if (eps->eigi[i] < 0) {
134:         eps->eigi[i] = -eps->eigi[i];
135:         eps->eigi[i+1] = -eps->eigi[i+1];
136:         if (!eps->evecsavailable) {
137:           /* the next correction only works with eigenvectors */
138:           (*eps->ops->computevectors)(eps);
139:         }
140:         VecScale(eps->V[i+1],-1.0);
141:       }
142:       i++;
143:     }
144:   }
145: #endif

147:   /* sort eigenvalues according to eps->which parameter */
148:   PetscFree(eps->perm);
149:   if (eps->nconv > 0) {
150:     PetscMalloc(sizeof(PetscInt)*eps->nconv, &eps->perm);
151:     EPSSortEigenvalues(eps->nconv, eps->eigr, eps->eigi, eps->which, eps->nconv, eps->perm);
152:   }

154:   PetscOptionsHasName(((PetscObject)eps)->prefix,"-eps_view",&flg);
155:   if (flg && !PetscPreLoadingOn) { EPSView(eps,PETSC_VIEWER_STDOUT_WORLD); }

157:   PetscOptionsHasName(((PetscObject)eps)->prefix,"-eps_plot_eigs",&flg);
158:   if (flg) {
159:     PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
160:                              PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
161:     PetscViewerDrawGetDraw(viewer,0,&draw);
162:     PetscDrawSPCreate(draw,1,&drawsp);
163:     for( i=0; i<eps->nconv; i++ ) {
164: #if defined(PETSC_USE_COMPLEX)
165:       re = PetscRealPart(eps->eigr[i]);
166:       im = PetscImaginaryPart(eps->eigi[i]);
167: #else
168:       re = eps->eigr[i];
169:       im = eps->eigi[i];
170: #endif
171:       PetscDrawSPAddPoint(drawsp,&re,&im);
172:     }
173:     PetscDrawSPDraw(drawsp);
174:     PetscDrawSPDestroy(drawsp);
175:     PetscViewerDestroy(viewer);
176:   }

178:   return(0);
179: }

183: /*@
184:    EPSGetIterationNumber - Gets the current iteration number. If the 
185:    call to EPSSolve() is complete, then it returns the number of iterations 
186:    carried out by the solution method.
187:  
188:    Not Collective

190:    Input Parameter:
191: .  eps - the eigensolver context

193:    Output Parameter:
194: .  its - number of iterations

196:    Level: intermediate

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

205: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
206: @*/
207: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
208: {
212:   *its = eps->its;
213:   return(0);
214: }

218: /*@
219:    EPSGetOperationCounters - Gets the total number of operator applications,
220:    inner product operations and linear iterations used by the ST object 
221:    during the last EPSSolve() call.

223:    Not Collective

225:    Input Parameter:
226: .  eps - EPS context

228:    Output Parameter:
229: +  ops  - number of operator applications
230: .  dots - number of inner product operations
231: -  lits - number of linear iterations

233:    Notes:
234:    When the eigensolver algorithm invokes STApply() then a linear system 
235:    must be solved (except in the case of standard eigenproblems and shift
236:    transformation). The number of iterations required in this solve is
237:    accumulated into a counter whose value is returned by this function.

239:    These counters are reset to zero at each successive call to EPSSolve().

241:    Level: intermediate

243: @*/
244: PetscErrorCode EPSGetOperationCounters(EPS eps,PetscInt* ops,PetscInt* dots,PetscInt* lits)
245: {

250:   STGetOperationCounters(eps->OP,ops,lits);
251:   if (dots) {
252:     IPGetOperationCounters(eps->ip,dots);
253:   }
254:   return(0);
255: }

259: /*@
260:    EPSGetConverged - Gets the number of converged eigenpairs.

262:    Not Collective

264:    Input Parameter:
265: .  eps - the eigensolver context
266:   
267:    Output Parameter:
268: .  nconv - number of converged eigenpairs 

270:    Note:
271:    This function should be called after EPSSolve() has finished.

273:    Level: beginner

275: .seealso: EPSSetDimensions()
276: @*/
277: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
278: {
282:   *nconv = eps->nconv;
283:   return(0);
284: }


289: /*@C
290:    EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was 
291:    stopped.

293:    Not Collective

295:    Input Parameter:
296: .  eps - the eigensolver context

298:    Output Parameter:
299: .  reason - negative value indicates diverged, positive value converged
300:    (see EPSConvergedReason)

302:    Possible values for reason:
303: +  EPS_CONVERGED_TOL - converged up to tolerance
304: .  EPS_DIVERGED_ITS - required more than its to reach convergence
305: .  EPS_DIVERGED_BREAKDOWN - generic breakdown in method
306: -  EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric

308:    Level: intermediate

310:    Notes: Can only be called after the call to EPSSolve() is complete.

312: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
313: @*/
314: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
315: {
319:   *reason = eps->reason;
320:   return(0);
321: }

325: /*@
326:    EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant 
327:    subspace.

329:    Not Collective

331:    Input Parameter:
332: .  eps - the eigensolver context
333:   
334:    Output Parameter:
335: .  v - an array of vectors

337:    Notes:
338:    This function should be called after EPSSolve() has finished.

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

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

350:    Level: intermediate

352: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetLeftInvariantSubspace()
353: @*/
354: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
355: {
357:   PetscInt       i;

363:   if (!eps->V) {
364:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
365:   }
366:   if (!eps->ishermitian && eps->evecsavailable) {
367:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetRightVector,EPSComputeRelativeError or EPSComputeResidualNorm");
368:   }
369:   for (i=0;i<eps->nconv;i++) {
370:     VecCopy(eps->V[i],v[i]);
371:   }
372:   return(0);
373: }

377: /*@
378:    EPSGetLeftInvariantSubspace - Gets an orthonormal basis of the computed left
379:    invariant subspace (only available in two-sided eigensolvers).

381:    Not Collective

383:    Input Parameter:
384: .  eps - the eigensolver context
385:   
386:    Output Parameter:
387: .  v - an array of vectors

389:    Notes:
390:    This function should be called after EPSSolve() has finished.

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

395:    The first k vectors returned in v span a left invariant subspace associated 
396:    with the first k computed eigenvalues (note that this is not true if the 
397:    k-th eigenvalue is complex and matrix A is real; in this case the first 
398:    k+1 vectors should be used). A left invariant subspace Y of A satisfies y'A 
399:    in Y for all y in Y (a similar definition applies for generalized 
400:    eigenproblems). 

402:    Level: intermediate

404: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
405: @*/
406: PetscErrorCode EPSGetLeftInvariantSubspace(EPS eps, Vec *v)
407: {
409:   PetscInt       i;

415:   if (!eps->W) {
416:     if (eps->solverclass!=EPS_TWO_SIDE) {
417:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
418:     } else {
419:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
420:     }
421:   }
422:   for (i=0;i<eps->nconv;i++) {
423:     VecCopy(eps->W[i],v[i]);
424:   }
425:   return(0);
426: }

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

434:    Not Collective

436:    Input Parameters:
437: +  eps - eigensolver context 
438: -  i   - index of the solution

440:    Output Parameters:
441: +  eigr - real part of eigenvalue
442: .  eigi - imaginary part of eigenvalue
443: .  Vr   - real part of eigenvector
444: -  Vi   - imaginary part of eigenvector

446:    Notes:
447:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is 
448:    configured with complex scalars the eigenvalue is stored 
449:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is 
450:    set to zero).

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

456:    Level: beginner

458: .seealso: EPSGetValue(), EPSGetRightVector(), EPSGetLeftVector(), EPSSolve(), 
459:           EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
460: @*/
461: PetscErrorCode EPSGetEigenpair(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
462: {

467:   if (!eps->eigr || !eps->eigi || !eps->V) {
468:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
469:   }
470:   if (i<0 || i>=eps->nconv) {
471:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
472:   }
473:   EPSGetValue(eps,i,eigr,eigi);
474:   EPSGetRightVector(eps,i,Vr,Vi);
475: 
476:   return(0);
477: }

481: /*@
482:    EPSGetValue - Gets the i-th eigenvalue as computed by EPSSolve(). 

484:    Not Collective

486:    Input Parameters:
487: +  eps - eigensolver context 
488: -  i   - index of the solution

490:    Output Parameters:
491: +  eigr - real part of eigenvalue
492: -  eigi - imaginary part of eigenvalue

494:    Notes:
495:    If the eigenvalue is real, then eigi is set to zero. If PETSc is 
496:    configured with complex scalars the eigenvalue is stored 
497:    directly in eigr (eigi is set to zero).

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

503:    Level: beginner

505: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), 
506:           EPSGetEigenpair()
507: @*/
508: PetscErrorCode EPSGetValue(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi)
509: {
510:   PetscInt       k;

514:   if (!eps->eigr || !eps->eigi) {
515:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
516:   }
517:   if (i<0 || i>=eps->nconv) {
518:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
519:   }

521:   if (!eps->perm) k = i;
522:   else k = eps->perm[i];
523: #ifdef PETSC_USE_COMPLEX
524:   if (eigr) *eigr = eps->eigr[k];
525:   if (eigi) *eigi = 0;
526: #else
527:   if (eigr) *eigr = eps->eigr[k];
528:   if (eigi) *eigi = eps->eigi[k];
529: #endif
530: 
531:   return(0);
532: }

536: /*@
537:    EPSGetRightVector - Gets the i-th right eigenvector as computed by EPSSolve(). 

539:    Not Collective

541:    Input Parameters:
542: +  eps - eigensolver context 
543: -  i   - index of the solution

545:    Output Parameters:
546: +  Vr   - real part of eigenvector
547: -  Vi   - imaginary part of eigenvector

549:    Notes:
550:    If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is 
551:    configured with complex scalars the eigenvector is stored 
552:    directly in Vr (Vi is set to zero).

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

558:    Level: beginner

560: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), 
561:           EPSGetEigenpair(), EPSGetLeftVector()
562: @*/
563: PetscErrorCode EPSGetRightVector(EPS eps, PetscInt i, Vec Vr, Vec Vi)
564: {
566:   PetscInt       k;

570:   if (!eps->V) {
571:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
572:   }
573:   if (i<0 || i>=eps->nconv) {
574:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
575:   }
576:   if (!eps->evecsavailable && (Vr || Vi) ) {
577:     (*eps->ops->computevectors)(eps);
578:   }

580:   if (!eps->perm) k = i;
581:   else k = eps->perm[i];
582: #ifdef PETSC_USE_COMPLEX
583:   if (Vr) { VecCopy(eps->V[k], Vr);  }
584:   if (Vi) { VecSet(Vi,0.0);  }
585: #else
586:   if (eps->eigi[k] > 0) { /* first value of conjugate pair */
587:     if (Vr) { VecCopy(eps->V[k], Vr);  }
588:     if (Vi) { VecCopy(eps->V[k+1], Vi);  }
589:   } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
590:     if (Vr) { VecCopy(eps->V[k-1], Vr);  }
591:     if (Vi) {
592:       VecCopy(eps->V[k], Vi);
593:       VecScale(Vi,-1.0);
594:     }
595:   } else { /* real eigenvalue */
596:     if (Vr) { VecCopy(eps->V[k], Vr);  }
597:     if (Vi) { VecSet(Vi,0.0);  }
598:   }
599: #endif
600: 
601:   return(0);
602: }

606: /*@
607:    EPSGetLeftVector - Gets the i-th left eigenvector as computed by EPSSolve() 
608:    (only available in two-sided eigensolvers). 

610:    Not Collective

612:    Input Parameters:
613: +  eps - eigensolver context 
614: -  i   - index of the solution

616:    Output Parameters:
617: +  Wr   - real part of eigenvector
618: -  Wi   - imaginary part of eigenvector

620:    Notes:
621:    If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is 
622:    configured with complex scalars the eigenvector is stored 
623:    directly in Wr (Wi is set to zero).

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

629:    Level: beginner

631: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), 
632:           EPSGetEigenpair(), EPSGetLeftVector()
633: @*/
634: PetscErrorCode EPSGetLeftVector(EPS eps, PetscInt i, Vec Wr, Vec Wi)
635: {
637:   PetscInt       k;

641:   if (!eps->W) {
642:     if (eps->solverclass!=EPS_TWO_SIDE) {
643:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
644:     } else {
645:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
646:     }
647:   }
648:   if (i<0 || i>=eps->nconv) {
649:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
650:   }
651:   if (!eps->evecsavailable && (Wr || Wi) ) {
652:     (*eps->ops->computevectors)(eps);
653:   }

655:   if (!eps->perm) k = i;
656:   else k = eps->perm[i];
657: #ifdef PETSC_USE_COMPLEX
658:   if (Wr) { VecCopy(eps->W[k], Wr);  }
659:   if (Wi) { VecSet(Wi,0.0);  }
660: #else
661:   if (eps->eigi[k] > 0) { /* first value of conjugate pair */
662:     if (Wr) { VecCopy(eps->W[k], Wr);  }
663:     if (Wi) { VecCopy(eps->W[k+1], Wi);  }
664:   } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
665:     if (Wr) { VecCopy(eps->W[k-1], Wr);  }
666:     if (Wi) {
667:       VecCopy(eps->W[k], Wi);
668:       VecScale(Wi,-1.0);
669:     }
670:   } else { /* real eigenvalue */
671:     if (Wr) { VecCopy(eps->W[k], Wr);  }
672:     if (Wi) { VecSet(Wi,0.0);  }
673:   }
674: #endif
675: 
676:   return(0);
677: }

681: /*@
682:    EPSGetErrorEstimate - Returns the error estimate associated to the i-th 
683:    computed eigenpair.

685:    Not Collective

687:    Input Parameter:
688: +  eps - eigensolver context 
689: -  i   - index of eigenpair

691:    Output Parameter:
692: .  errest - the error estimate

694:    Notes:
695:    This is the error estimate used internally by the eigensolver. The actual
696:    error bound can be computed with EPSComputeRelativeError(). See also the user's
697:    manual for details.

699:    Level: advanced

701: .seealso: EPSComputeRelativeError()
702: @*/
703: PetscErrorCode EPSGetErrorEstimate(EPS eps, PetscInt i, PetscReal *errest)
704: {
707:   if (!eps->eigr || !eps->eigi) {
708:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
709:   }
710:   if (i<0 || i>=eps->nconv) {
711:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
712:   }
713:   if (eps->perm) i = eps->perm[i];
714:   if (errest) *errest = eps->errest[i];
715:   return(0);
716: }

720: /*@
721:    EPSGetErrorEstimateLeft - Returns the left error estimate associated to the i-th 
722:    computed eigenpair (only available in two-sided eigensolvers).

724:    Not Collective

726:    Input Parameter:
727: +  eps - eigensolver context 
728: -  i   - index of eigenpair

730:    Output Parameter:
731: .  errest - the left error estimate

733:    Notes:
734:    This is the error estimate used internally by the eigensolver. The actual
735:    error bound can be computed with EPSComputeRelativeErrorLeft(). See also the user's
736:    manual for details.

738:    Level: advanced

740: .seealso: EPSComputeRelativeErrorLeft()
741: @*/
742: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps, PetscInt i, PetscReal *errest)
743: {
746:   if (!eps->eigr || !eps->eigi) {
747:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
748:   }
749:   if (eps->solverclass!=EPS_TWO_SIDE) {
750:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Only available for two-sided solvers");
751:   }
752:   if (i<0 || i>=eps->nconv) {
753:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
754:   }
755:   if (eps->perm) i = eps->perm[i];
756:   if (errest) *errest = eps->errest_left[i];
757:   return(0);
758: }

762: /*@
763:    EPSComputeResidualNorm - Computes the norm of the residual vector associated with 
764:    the i-th computed eigenpair.

766:    Collective on EPS

768:    Input Parameter:
769: .  eps - the eigensolver context
770: .  i   - the solution index

772:    Output Parameter:
773: .  norm - the residual norm, computed as ||Ax-kBx||_2 where k is the 
774:    eigenvalue and x is the eigenvector. 
775:    If k=0 then the residual norm is computed as ||Ax||_2.

777:    Notes:
778:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
779:    Eigenpairs are indexed according to the ordering criterion established 
780:    with EPSSetWhichEigenpairs().

782:    Level: beginner

784: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
785: @*/
786: PetscErrorCode EPSComputeResidualNorm(EPS eps, PetscInt i, PetscReal *norm)
787: {
789:   Vec            u, v, w, xr, xi;
790:   Mat            A, B;
791:   PetscScalar    kr, ki;
792: #ifndef PETSC_USE_COMPLEX
793:   PetscReal      ni, nr;
794: #endif
795: 
798:   STGetOperators(eps->OP,&A,&B);
799:   VecDuplicate(eps->vec_initial,&u);
800:   VecDuplicate(eps->vec_initial,&v);
801:   VecDuplicate(eps->vec_initial,&w);
802:   VecDuplicate(eps->vec_initial,&xr);
803:   VecDuplicate(eps->vec_initial,&xi);
804:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);

806: #ifndef PETSC_USE_COMPLEX
807:   if (ki == 0 ||
808:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
809: #endif
810:     MatMult( A, xr, u );  /* u=A*x */
811:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
812:       if (eps->isgeneralized) { MatMult( B, xr, w );  }
813:       else { VecCopy( xr, w );  } /* w=B*x */
814:       VecAXPY( u, -kr, w );  /* u=A*x-k*B*x */
815:     }
816:     VecNorm( u, NORM_2, norm);
817: #ifndef PETSC_USE_COMPLEX
818:   } else {
819:     MatMult( A, xr, u );  /* u=A*xr */
820:     if (eps->isgeneralized) { MatMult( B, xr, v );  }
821:     else { VecCopy( xr, v );  } /* v=B*xr */
822:     VecAXPY( u, -kr, v );  /* u=A*xr-kr*B*xr */
823:     if (eps->isgeneralized) { MatMult( B, xi, w );  }
824:     else { VecCopy( xi, w );  } /* w=B*xi */
825:     VecAXPY( u, ki, w );  /* u=A*xr-kr*B*xr+ki*B*xi */
826:     VecNorm( u, NORM_2, &nr );
827:     MatMult( A, xi, u );  /* u=A*xi */
828:     VecAXPY( u, -kr, w );  /* u=A*xi-kr*B*xi */
829:     VecAXPY( u, -ki, v );  /* u=A*xi-kr*B*xi-ki*B*xr */
830:     VecNorm( u, NORM_2, &ni );
831:     *norm = SlepcAbsEigenvalue( nr, ni );
832:   }
833: #endif

835:   VecDestroy(w);
836:   VecDestroy(v);
837:   VecDestroy(u);
838:   VecDestroy(xr);
839:   VecDestroy(xi);
840:   return(0);
841: }

845: /*@
846:    EPSComputeResidualNormLeft - Computes the norm of the residual vector associated with 
847:    the i-th computed left eigenvector (only available in two-sided eigensolvers).

849:    Collective on EPS

851:    Input Parameter:
852: .  eps - the eigensolver context
853: .  i   - the solution index

855:    Output Parameter:
856: .  norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the 
857:    eigenvalue and y is the left eigenvector. 
858:    If k=0 then the residual norm is computed as ||y'A||_2.

860:    Notes:
861:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
862:    Eigenpairs are indexed according to the ordering criterion established 
863:    with EPSSetWhichEigenpairs().

865:    Level: beginner

867: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
868: @*/
869: PetscErrorCode EPSComputeResidualNormLeft(EPS eps, PetscInt i, PetscReal *norm)
870: {
872:   Vec            u, v, w, xr, xi;
873:   Mat            A, B;
874:   PetscScalar    kr, ki;
875: #ifndef PETSC_USE_COMPLEX
876:   PetscReal      ni, nr;
877: #endif
878: 
881:   STGetOperators(eps->OP,&A,&B);
882:   VecDuplicate(eps->vec_initial_left,&u);
883:   VecDuplicate(eps->vec_initial_left,&v);
884:   VecDuplicate(eps->vec_initial_left,&w);
885:   VecDuplicate(eps->vec_initial_left,&xr);
886:   VecDuplicate(eps->vec_initial_left,&xi);
887:   EPSGetValue(eps,i,&kr,&ki);
888:   EPSGetLeftVector(eps,i,xr,xi);

890: #ifndef PETSC_USE_COMPLEX
891:   if (ki == 0 ||
892:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
893: #endif
894:     MatMultTranspose( A, xr, u );  /* u=A'*x */
895:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
896:       if (eps->isgeneralized) { MatMultTranspose( B, xr, w );  }
897:       else { VecCopy( xr, w );  } /* w=B'*x */
898:       VecAXPY( u, -kr, w);  /* u=A'*x-k*B'*x */
899:     }
900:     VecNorm( u, NORM_2, norm);
901: #ifndef PETSC_USE_COMPLEX
902:   } else {
903:     MatMultTranspose( A, xr, u );  /* u=A'*xr */
904:     if (eps->isgeneralized) { MatMultTranspose( B, xr, v );  }
905:     else { VecCopy( xr, v );  } /* v=B'*xr */
906:     VecAXPY( u, -kr, v );  /* u=A'*xr-kr*B'*xr */
907:     if (eps->isgeneralized) { MatMultTranspose( B, xi, w );  }
908:     else { VecCopy( xi, w );  } /* w=B'*xi */
909:     VecAXPY( u, ki, w );  /* u=A'*xr-kr*B'*xr+ki*B'*xi */
910:     VecNorm( u, NORM_2, &nr );
911:     MatMultTranspose( A, xi, u );  /* u=A'*xi */
912:     VecAXPY( u, -kr, w );  /* u=A'*xi-kr*B'*xi */
913:     VecAXPY( u, -ki, v );  /* u=A'*xi-kr*B'*xi-ki*B'*xr */
914:     VecNorm( u, NORM_2, &ni );
915:     *norm = SlepcAbsEigenvalue( nr, ni );
916:   }
917: #endif

919:   VecDestroy(w);
920:   VecDestroy(v);
921:   VecDestroy(u);
922:   VecDestroy(xr);
923:   VecDestroy(xi);
924:   return(0);
925: }

929: /*@
930:    EPSComputeRelativeError - Computes the relative error bound associated 
931:    with the i-th computed eigenpair.

933:    Collective on EPS

935:    Input Parameter:
936: .  eps - the eigensolver context
937: .  i   - the solution index

939:    Output Parameter:
940: .  error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where 
941:    k is the eigenvalue and x is the eigenvector. 
942:    If k=0 the relative error is computed as ||Ax||_2/||x||_2.

944:    Level: beginner

946: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
947: @*/
948: PetscErrorCode EPSComputeRelativeError(EPS eps, PetscInt i, PetscReal *error)
949: {
951:   Vec            xr, xi;
952:   PetscScalar    kr, ki;
953:   PetscReal      norm, er;
954: #ifndef PETSC_USE_COMPLEX
955:   Vec            u;
956:   PetscReal      ei;
957: #endif
958: 
961:   EPSComputeResidualNorm(eps,i,&norm);
962:   VecDuplicate(eps->vec_initial,&xr);
963:   VecDuplicate(eps->vec_initial,&xi);
964:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);

966: #ifndef PETSC_USE_COMPLEX
967:   if (ki == 0 ||
968:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
969: #endif
970:     VecNorm(xr, NORM_2, &er);
971:     if (PetscAbsScalar(kr) > norm) {
972:       *error =  norm / (PetscAbsScalar(kr) * er);
973:     } else {
974:       *error = norm / er;
975:     }
976: #ifndef PETSC_USE_COMPLEX
977:   } else {
978:     if (SlepcAbsEigenvalue(kr,ki) > norm) {
979:       VecDuplicate(xi, &u);
980:       VecCopy(xi, u);
981:       VecAXPBY(u, kr, -ki, xr);
982:       VecNorm(u, NORM_2, &er);
983:       VecAXPBY(xi, kr, ki, xr);
984:       VecNorm(xi, NORM_2, &ei);
985:       VecDestroy(u);
986:     } else {
987:       VecDot(xr, xr, &er);
988:       VecDot(xi, xi, &ei);
989:     }
990:     *error = norm / SlepcAbsEigenvalue(er, ei);
991:   }
992: #endif    
993: 
994:   VecDestroy(xr);
995:   VecDestroy(xi);
996:   return(0);
997: }

1001: /*@
1002:    EPSComputeRelativeErrorLeft - Computes the relative error bound associated 
1003:    with the i-th computed eigenvalue and left eigenvector (only available in
1004:    two-sided eigensolvers).

1006:    Collective on EPS

1008:    Input Parameter:
1009: .  eps - the eigensolver context
1010: .  i   - the solution index

1012:    Output Parameter:
1013: .  error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where 
1014:    k is the eigenvalue and y is the left eigenvector. 
1015:    If k=0 the relative error is computed as ||y'A||_2/||y||_2.

1017:    Level: beginner

1019: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
1020: @*/
1021: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps, PetscInt i, PetscReal *error)
1022: {
1024:   Vec            xr, xi;
1025:   PetscScalar    kr, ki;
1026:   PetscReal      norm, er;
1027: #ifndef PETSC_USE_COMPLEX
1028:   Vec            u;
1029:   PetscReal      ei;
1030: #endif
1031: 
1034:   EPSComputeResidualNormLeft(eps,i,&norm);
1035:   VecDuplicate(eps->vec_initial_left,&xr);
1036:   VecDuplicate(eps->vec_initial_left,&xi);
1037:   EPSGetValue(eps,i,&kr,&ki);
1038:   EPSGetLeftVector(eps,i,xr,xi);

1040: #ifndef PETSC_USE_COMPLEX
1041:   if (ki == 0 ||
1042:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1043: #endif
1044:     VecNorm(xr, NORM_2, &er);
1045:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1046:       *error =  norm / (PetscAbsScalar(kr) * er);
1047:     } else {
1048:       *error = norm / er;
1049:     }
1050: #ifndef PETSC_USE_COMPLEX
1051:   } else {
1052:     VecDuplicate(xi, &u);
1053:     VecCopy(xi, u);
1054:     VecAXPBY(u, kr, -ki, xr);
1055:     VecNorm(u, NORM_2, &er);
1056:     VecAXPBY(xi, kr, ki, xr);
1057:     VecNorm(xi, NORM_2, &ei);
1058:     VecDestroy(u);
1059:     *error = norm / SlepcAbsEigenvalue(er, ei);
1060:   }
1061: #endif    
1062: 
1063:   VecDestroy(xr);
1064:   VecDestroy(xi);
1065:   return(0);
1066: }

1068: #define SWAP(a,b,t) {t=a;a=b;b=t;}

1072: /*@
1073:    EPSSortEigenvalues - Sorts a list of eigenvalues according to a certain
1074:    criterion.

1076:    Not Collective

1078:    Input Parameters:
1079: +  n     - number of eigenvalue in the list
1080: .  eig   - pointer to the array containing the eigenvalues
1081: .  eigi  - imaginary part of the eigenvalues (only when using real numbers)
1082: .  which - sorting criterion
1083: -  nev   - number of wanted eigenvalues

1085:    Output Parameter:
1086: .  permout - resulting permutation

1088:    Notes:
1089:    The result is a list of indices in the original eigenvalue array 
1090:    corresponding to the first nev eigenvalues sorted in the specified
1091:    criterion

1093:    Level: developer

1095: .seealso: EPSSortEigenvaluesReal(), EPSSetWhichEigenpairs()
1096: @*/
1097: PetscErrorCode EPSSortEigenvalues(PetscInt n,PetscScalar *eig,PetscScalar *eigi,EPSWhich which,PetscInt nev,PetscInt *permout)
1098: {
1100:   PetscInt       i;
1101:   PetscInt       *perm;
1102:   PetscReal      *values;

1105:   PetscMalloc(n*sizeof(PetscInt),&perm);
1106:   PetscMalloc(n*sizeof(PetscReal),&values);
1107:   for (i=0; i<n; i++) { perm[i] = i;}

1109:   switch(which) {
1110:     case EPS_LARGEST_MAGNITUDE:
1111:     case EPS_SMALLEST_MAGNITUDE:
1112:       for (i=0; i<n; i++) { values[i] = SlepcAbsEigenvalue(eig[i],eigi[i]); }
1113:       break;
1114:     case EPS_LARGEST_REAL:
1115:     case EPS_SMALLEST_REAL:
1116:       for (i=0; i<n; i++) { values[i] = PetscRealPart(eig[i]); }
1117:       break;
1118:     case EPS_LARGEST_IMAGINARY:
1119:     case EPS_SMALLEST_IMAGINARY:
1120: #if defined(PETSC_USE_COMPLEX)
1121:       for (i=0; i<n; i++) { values[i] = PetscImaginaryPart(eig[i]); }
1122: #else
1123:       for (i=0; i<n; i++) { values[i] = PetscAbsReal(eigi[i]); }
1124: #endif
1125:       break;
1126:     default: SETERRQ(1,"Wrong value of which");
1127:   }

1129:   PetscSortRealWithPermutation(n,values,perm);

1131:   switch(which) {
1132:     case EPS_LARGEST_MAGNITUDE:
1133:     case EPS_LARGEST_REAL:
1134:     case EPS_LARGEST_IMAGINARY:
1135:       for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
1136:       break;
1137:     case EPS_SMALLEST_MAGNITUDE:
1138:     case EPS_SMALLEST_REAL:
1139:     case EPS_SMALLEST_IMAGINARY:
1140:       for (i=0; i<nev; i++) { permout[i] = perm[i]; }
1141:       break;
1142:     default: SETERRQ(1,"Wrong value of which");
1143:   }

1145: #if !defined(PETSC_USE_COMPLEX)
1146:   for (i=0; i<nev-1; i++) {
1147:     if (eigi[permout[i]] != 0.0) {
1148:       if (eig[permout[i]] == eig[permout[i+1]] &&
1149:           eigi[permout[i]] == -eigi[permout[i+1]] &&
1150:           eigi[permout[i]] < 0.0) {
1151:         PetscInt tmp;
1152:         SWAP(permout[i], permout[i+1], tmp);
1153:       }
1154:       i++;
1155:     }
1156:   }
1157: #endif

1159:   PetscFree(values);
1160:   PetscFree(perm);
1161:   return(0);
1162: }

1166: /*@
1167:    EPSSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
1168:    criterion (version for real eigenvalues only).

1170:    Not Collective

1172:    Input Parameters:
1173: +  n     - number of eigenvalue in the list
1174: .  eig   - pointer to the array containing the eigenvalues (real)
1175: .  which - sorting criterion
1176: -  nev   - number of wanted eigenvalues

1178:    Output Parameter:
1179: .  permout - resulting permutation

1181:    Workspace:
1182: .  work - workspace for storing n real values and n integer values

1184:    Notes:
1185:    The result is a list of indices in the original eigenvalue array 
1186:    corresponding to the first nev eigenvalues sorted in the specified
1187:    criterion

1189:    Level: developer

1191: .seealso: EPSSortEigenvalues(), EPSSetWhichEigenpairs()
1192: @*/
1193: PetscErrorCode EPSSortEigenvaluesReal(PetscInt n,PetscReal *eig,EPSWhich which,PetscInt nev,PetscInt *permout,PetscReal *work)
1194: {
1196:   PetscInt            i;
1197:   PetscReal      *values = work;
1198:   PetscInt       *perm = (PetscInt*)(work+n);

1201:   for (i=0; i<n; i++) { perm[i] = i;}

1203:   switch(which) {
1204:     case EPS_LARGEST_MAGNITUDE:
1205:     case EPS_SMALLEST_MAGNITUDE:
1206:       for (i=0; i<n; i++) { values[i] = PetscAbsReal(eig[i]); }
1207:       break;
1208:     case EPS_LARGEST_REAL:
1209:     case EPS_SMALLEST_REAL:
1210:       for (i=0; i<n; i++) { values[i] = eig[i]; }
1211:       break;
1212:     default: SETERRQ(1,"Wrong value of which");
1213:   }

1215:   PetscSortRealWithPermutation(n,values,perm);

1217:   switch(which) {
1218:     case EPS_LARGEST_MAGNITUDE:
1219:     case EPS_LARGEST_REAL:
1220:       for (i=0; i<nev; i++) { permout[i] = perm[n-1-i]; }
1221:       break;
1222:     case EPS_SMALLEST_MAGNITUDE:
1223:     case EPS_SMALLEST_REAL:
1224:       for (i=0; i<nev; i++) { permout[i] = perm[i]; }
1225:       break;
1226:     default: SETERRQ(1,"Wrong value of which");
1227:   }
1228:   return(0);
1229: }

1233: /*@
1234:    EPSGetStartVector - Gets a vector to be used as the starting vector
1235:    in an Arnoldi or Lanczos reduction.

1237:    Collective on EPS and Vec

1239:    Input Parameters:
1240: +  eps - the eigensolver context
1241: -  i   - index of the Arnoldi/Lanczos step

1243:    Output Parameters:
1244: +  vec - the start vector
1245: -  breakdown - flag indicating that a breakdown has occurred

1247:    Notes:
1248:    The start vector is computed from another vector: for the first step (i=0),
1249:    the initial vector is used (see EPSGetInitialVector()); otherwise a random
1250:    vector is created. Then this vector is forced to be in the range of OP (only
1251:    for generalized definite problems) and orthonormalized with respect to all
1252:    V-vectors up to i-1.

1254:    The flag breakdown is set to true if either i=0 and the vector belongs to the
1255:    deflation space, or i>0 and the vector is linearly dependent with respect
1256:    to the V-vectors.

1258:    The caller must pass a vector already allocated with dimensions conforming
1259:    to the initial vector. This vector is overwritten.

1261:    Level: developer

1263: .seealso: EPSGetInitialVector()

1265: @*/
1266: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,Vec vec,PetscTruth *breakdown)
1267: {
1269:   PetscReal      norm;
1270:   PetscTruth     lindep;
1271:   Vec            w;
1272: 

1277:   /* For the first step, use the initial vector, otherwise a random one */
1278:   if (i==0) {
1279:     w = eps->vec_initial;
1280:   } else {
1281:     VecDuplicate(eps->vec_initial,&w);
1282:     SlepcVecSetRandom(w);
1283:   }

1285:   /* Force the vector to be in the range of OP for definite generalized problems */
1286:   if (eps->ispositive) {
1287:     STApply(eps->OP,w,vec);
1288:   } else {
1289:     VecCopy(w,vec);
1290:   }

1292:   /* Orthonormalize the vector with respect to previous vectors */
1293:   IPOrthogonalize(eps->ip,i+eps->nds,PETSC_NULL,eps->DSV,vec,PETSC_NULL,&norm,&lindep,PETSC_NULL,PETSC_NULL);
1294:   if (breakdown) *breakdown = lindep;
1295:   else if (lindep || norm == 0.0) {
1296:     if (i==0) { SETERRQ(1,"Initial vector is zero or belongs to the deflation space"); }
1297:     else { SETERRQ(1,"Unable to generate more start vectors"); }
1298:   }
1299: 
1300:   VecScale(vec,1.0/norm);

1302:   if (i!=0) {
1303:     VecDestroy(w);
1304:   }

1306:   return(0);
1307: }
1310: /*@
1311:    EPSGetLeftStartVector - Gets a vector to be used as the starting vector
1312:    in the left recurrence of a two-sided eigensolver.

1314:    Collective on EPS and Vec

1316:    Input Parameters:
1317: +  eps - the eigensolver context
1318: -  i   - index of the Arnoldi/Lanczos step

1320:    Output Parameter:
1321: .  vec - the start vector

1323:    Notes:
1324:    The start vector is computed from another vector: for the first step (i=0),
1325:    the left initial vector is used (see EPSGetLeftInitialVector()); otherwise 
1326:    a random vector is created. Then this vector is forced to be in the range 
1327:    of OP' and orthonormalized with respect to all W-vectors up to i-1.

1329:    The caller must pass a vector already allocated with dimensions conforming
1330:    to the left initial vector. This vector is overwritten.

1332:    Level: developer

1334: .seealso: EPSGetLeftInitialVector()

1336: @*/
1337: PetscErrorCode EPSGetLeftStartVector(EPS eps,PetscInt i,Vec vec)
1338: {
1340:   PetscTruth     breakdown;
1341:   PetscReal      norm;
1342:   Vec            w;
1343: 

1348:   /* For the first step, use the initial vector, otherwise a random one */
1349:   if (i==0) {
1350:     w = eps->vec_initial_left;
1351:   }
1352:   else {
1353:     VecDuplicate(eps->vec_initial_left,&w);
1354:     SlepcVecSetRandom(w);
1355:   }

1357:   /* Force the vector to be in the range of OP */
1358:   STApplyTranspose(eps->OP,w,vec);

1360:   /* Orthonormalize the vector with respect to previous vectors */
1361:   IPOrthogonalize(eps->ip,i,PETSC_NULL,eps->W,vec,PETSC_NULL,&norm,&breakdown,PETSC_NULL,PETSC_NULL);
1362:   if (breakdown) {
1363:     if (i==0) { SETERRQ(1,"Left initial vector is zero"); }
1364:     else { SETERRQ(1,"Unable to generate more left start vectors"); }
1365:   }
1366:   VecScale(vec,1/norm);

1368:   if (i!=0) {
1369:     VecDestroy(w);
1370:   }

1372:   return(0);
1373: }