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-2010, 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

 26: typedef struct {
 27:   /* old values of eps */
 28:   EPSWhich
 29:     old_which;
 30:   PetscErrorCode
 31:     (*old_which_func)(EPS,PetscScalar,PetscScalar,PetscScalar,PetscScalar,
 32:                       PetscInt*,void*);
 33:   void
 34:     *old_which_ctx;
 35: } EPSSortForSTData;


 38: PetscErrorCode EPSSortForSTFunc(EPS eps, PetscScalar ar, PetscScalar ai,
 39:                                 PetscScalar br, PetscScalar bi, PetscInt *r,
 40:                                 void *ctx)
 41: {
 42:   EPSSortForSTData *data = (EPSSortForSTData*)ctx;
 43:   PetscErrorCode  ierr;


 47:   /* Back-transform the harmonic values */
 48:   STBackTransform(eps->OP,1,&ar,&ai);
 49:   STBackTransform(eps->OP,1,&br,&bi);

 51:   /* Compare values using the user options for the eigenpairs selection */
 52:   eps->which = data->old_which;
 53:   eps->which_func = data->old_which_func;
 54:   eps->which_ctx = data->old_which_ctx;
 55:   EPSCompareEigenvalues(eps, ar, ai, br, bi, r);

 57:   /* Restore the eps values */
 58:   eps->which = EPS_WHICH_USER;
 59:   eps->which_func = EPSSortForSTFunc;
 60:   eps->which_ctx = data;

 62:   return(0);
 63: }


 68: /*@
 69:    EPSSolve - Solves the eigensystem.

 71:    Collective on EPS

 73:    Input Parameter:
 74: .  eps - eigensolver context obtained from EPSCreate()

 76:    Options Database:
 77: +   -eps_view - print information about the solver used
 78: .   -eps_view_binary - save the matrices to the default binary file
 79: -   -eps_plot_eigs - plot computed eigenvalues

 81:    Level: beginner

 83: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances() 
 84: @*/
 85: PetscErrorCode EPSSolve(EPS eps)
 86: {
 88:   PetscInt       i;
 89:   PetscReal      re,im;
 90:   PetscScalar    dot;
 91:   PetscTruth     flg,isfold;
 92:   PetscViewer    viewer;
 93:   PetscDraw      draw;
 94:   PetscDrawSP    drawsp;
 95:   STMatMode      matmode;
 96:   char           filename[PETSC_MAX_PATH_LEN];
 97:   EPSSortForSTData data;
 98:   Mat            A,B;
 99:   KSP            ksp;
100:   Vec            w,x;


105:   flg = PETSC_FALSE;
106:   PetscOptionsGetTruth(((PetscObject)eps)->prefix,"-eps_view_binary",&flg,PETSC_NULL);
107:   if (flg) {
108:     STGetOperators(eps->OP,&A,&B);
109:     MatView(A,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
110:     if (B) MatView(B,PETSC_VIEWER_BINARY_(((PetscObject)eps)->comm));
111:   }

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

116:   /* call setup */
117:   if (!eps->setupcalled){ EPSSetUp(eps); }
118:   STResetOperationCounters(eps->OP);
119:   IPResetOperationCounters(eps->ip);
120:   eps->evecsavailable = PETSC_FALSE;
121:   eps->nconv = 0;
122:   eps->its = 0;
123:   for (i=0;i<eps->ncv;i++) eps->eigr[i]=eps->eigi[i]=eps->errest[i]=0.0;
124:   EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,eps->ncv);

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

128:   /* temporarily change which */
129:   data.old_which = eps->which;
130:   data.old_which_func = eps->which_func;
131:   data.old_which_ctx = eps->which_ctx;
132:   eps->which = EPS_WHICH_USER;
133:   eps->which_func = EPSSortForSTFunc;
134:   eps->which_ctx = &data;

136:   /* call solver */
137:   (*eps->ops->solve)(eps);

139:   /* restore which */
140:   eps->which = data.old_which;
141:   eps->which_func = data.old_which_func;
142:   eps->which_ctx = data.old_which_ctx;

144:   STGetMatMode(eps->OP,&matmode);
145:   if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
146:     /* Purify eigenvectors before reverting operator */
147:     (*eps->ops->computevectors)(eps);
148:   }
149:   STPostSolve(eps->OP);
150:   PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);

152:   if (!eps->reason) {
153:     SETERRQ(1,"Internal error, solver returned without setting converged reason");
154:   }

156:   /* Map eigenvalues back to the original problem, necessary in some 
157:   * spectral transformations */
158:   if (eps->ops->backtransform) {
159:     (*eps->ops->backtransform)(eps);
160:   }

162:   /* Adjust left eigenvectors in generalized problems: y = B^T y */
163:   if (eps->isgeneralized && eps->leftvecs) {
164:     STGetOperators(eps->OP,PETSC_NULL,&B);
165:     KSPCreate(((PetscObject)eps)->comm,&ksp);
166:     KSPSetOperators(ksp,B,B,DIFFERENT_NONZERO_PATTERN);
167:     KSPSetFromOptions(ksp);
168:     MatGetVecs(B,PETSC_NULL,&w);
169:     for (i=0;i<eps->nconv;i++) {
170:       VecCopy(eps->W[i],w);
171:       KSPSolveTranspose(ksp,w,eps->W[i]);
172:     }
173:     KSPDestroy(ksp);
174:     VecDestroy(w);
175:   }

177: #ifndef PETSC_USE_COMPLEX
178:   /* reorder conjugate eigenvalues (positive imaginary first) */
179:   for (i=0; i<eps->nconv-1; i++) {
180:     if (eps->eigi[i] != 0) {
181:       if (eps->eigi[i] < 0) {
182:         eps->eigi[i] = -eps->eigi[i];
183:         eps->eigi[i+1] = -eps->eigi[i+1];
184:         if (!eps->evecsavailable) {
185:           /* the next correction only works with eigenvectors */
186:           (*eps->ops->computevectors)(eps);
187:         }
188:         VecScale(eps->V[i+1],-1.0);
189:       }
190:       i++;
191:     }
192:   }
193: #endif

195:   /* quick and dirty solution for FOLD: recompute eigenvalues as Rayleigh quotients */
196:   PetscTypeCompare((PetscObject)eps->OP,STFOLD,&isfold);
197:   if (isfold) {
198:     STGetOperators(eps->OP,&A,&B);
199:     MatGetVecs(A,&w,PETSC_NULL);
200:     if (!eps->evecsavailable) { (*eps->ops->computevectors)(eps); }
201:     for (i=0;i<eps->nconv;i++) {
202:       x = eps->V[i];
203:       MatMult(A,x,w);
204:       VecDot(w,x,&eps->eigr[i]);
205:       if (eps->isgeneralized) {
206:         MatMult(B,x,w);
207:         VecDot(w,x,&dot);
208:         eps->eigr[i] /= dot;
209:       }
210:     }
211:     VecDestroy(w);
212:   }

214:   /* sort eigenvalues according to eps->which parameter */
215:   PetscFree(eps->perm);
216:   if (eps->nconv > 0) {
217:     PetscMalloc(sizeof(PetscInt)*eps->nconv, &eps->perm);
218:     EPSSortEigenvalues(eps, eps->nconv, eps->eigr, eps->eigi, eps->perm);
219:   }

221:   PetscOptionsGetString(((PetscObject)eps)->prefix,"-eps_view",filename,PETSC_MAX_PATH_LEN,&flg);
222:   if (flg && !PetscPreLoadingOn) {
223:     PetscViewerASCIIOpen(((PetscObject)eps)->comm,filename,&viewer);
224:     EPSView(eps,viewer);
225:     PetscViewerDestroy(viewer);
226:   }

228:   flg = PETSC_FALSE;
229:   PetscOptionsGetTruth(((PetscObject)eps)->prefix,"-eps_plot_eigs",&flg,PETSC_NULL);
230:   if (flg) {
231:     PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
232:                              PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
233:     PetscViewerDrawGetDraw(viewer,0,&draw);
234:     PetscDrawSPCreate(draw,1,&drawsp);
235:     for( i=0; i<eps->nconv; i++ ) {
236: #if defined(PETSC_USE_COMPLEX)
237:       re = PetscRealPart(eps->eigr[i]);
238:       im = PetscImaginaryPart(eps->eigi[i]);
239: #else
240:       re = eps->eigr[i];
241:       im = eps->eigi[i];
242: #endif
243:       PetscDrawSPAddPoint(drawsp,&re,&im);
244:     }
245:     PetscDrawSPDraw(drawsp);
246:     PetscDrawSPDestroy(drawsp);
247:     PetscViewerDestroy(viewer);
248:   }

250:   /* Remove the initial subspaces */
251:   eps->nini = 0;
252:   eps->ninil = 0;

254:   return(0);
255: }

259: /*@
260:    EPSGetIterationNumber - Gets the current iteration number. If the 
261:    call to EPSSolve() is complete, then it returns the number of iterations 
262:    carried out by the solution method.
263:  
264:    Not Collective

266:    Input Parameter:
267: .  eps - the eigensolver context

269:    Output Parameter:
270: .  its - number of iterations

272:    Level: intermediate

274:    Note:
275:    During the i-th iteration this call returns i-1. If EPSSolve() is 
276:    complete, then parameter "its" contains either the iteration number at
277:    which convergence was successfully reached, or failure was detected.  
278:    Call EPSGetConvergedReason() to determine if the solver converged or 
279:    failed and why.

281: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
282: @*/
283: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
284: {
288:   *its = eps->its;
289:   return(0);
290: }

294: /*@
295:    EPSGetOperationCounters - Gets the total number of operator applications,
296:    inner product operations and linear iterations used by the ST object 
297:    during the last EPSSolve() call.

299:    Not Collective

301:    Input Parameter:
302: .  eps - EPS context

304:    Output Parameter:
305: +  ops  - number of operator applications
306: .  dots - number of inner product operations
307: -  lits - number of linear iterations

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

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

317:    Level: intermediate

319: @*/
320: PetscErrorCode EPSGetOperationCounters(EPS eps,PetscInt* ops,PetscInt* dots,PetscInt* lits)
321: {

326:   STGetOperationCounters(eps->OP,ops,lits);
327:   if (dots) {
328:     IPGetOperationCounters(eps->ip,dots);
329:   }
330:   return(0);
331: }

335: /*@
336:    EPSGetConverged - Gets the number of converged eigenpairs.

338:    Not Collective

340:    Input Parameter:
341: .  eps - the eigensolver context
342:   
343:    Output Parameter:
344: .  nconv - number of converged eigenpairs 

346:    Note:
347:    This function should be called after EPSSolve() has finished.

349:    Level: beginner

351: .seealso: EPSSetDimensions(), EPSSolve()
352: @*/
353: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
354: {
358:   *nconv = eps->nconv;
359:   return(0);
360: }


365: /*@C
366:    EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was 
367:    stopped.

369:    Not Collective

371:    Input Parameter:
372: .  eps - the eigensolver context

374:    Output Parameter:
375: .  reason - negative value indicates diverged, positive value converged

377:    Possible values for reason:
378: +  EPS_CONVERGED_TOL - converged up to tolerance
379: .  EPS_DIVERGED_ITS - required more than its to reach convergence
380: .  EPS_DIVERGED_BREAKDOWN - generic breakdown in method
381: -  EPS_DIVERGED_NONSYMMETRIC - The operator is nonsymmetric

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

386:    Level: intermediate

388: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
389: @*/
390: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
391: {
395:   *reason = eps->reason;
396:   return(0);
397: }

401: /*@
402:    EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant 
403:    subspace.

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

407:    Input Parameter:
408: .  eps - the eigensolver context
409:   
410:    Output Parameter:
411: .  v - an array of vectors

413:    Notes:
414:    This function should be called after EPSSolve() has finished.

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

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

426:    Level: intermediate

428: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspaceLeft()
429: @*/
430: PetscErrorCode EPSGetInvariantSubspace(EPS eps, Vec *v)
431: {
433:   PetscInt       i;

439:   if (!eps->V) {
440:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
441:   }
442:   if (!eps->ishermitian && eps->evecsavailable) {
443:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector,EPSComputeRelativeError or EPSComputeResidualNorm");
444:   }
445:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
446:     for (i=0;i<eps->nconv;i++) {
447:       VecPointwiseDivide(v[i],eps->V[i],eps->D);
448:       VecNormalize(v[i],PETSC_NULL);
449:     }
450:   }
451:   else {
452:     for (i=0;i<eps->nconv;i++) {
453:       VecCopy(eps->V[i],v[i]);
454:     }
455:   }
456:   return(0);
457: }

461: /*@
462:    EPSGetInvariantSubspaceLeft - Gets an orthonormal basis of the computed left
463:    invariant subspace (only available in two-sided eigensolvers).

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

467:    Input Parameter:
468: .  eps - the eigensolver context
469:   
470:    Output Parameter:
471: .  v - an array of vectors

473:    Notes:
474:    This function should be called after EPSSolve() has finished.

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

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

486:    Level: intermediate

488: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve(), EPSGetInvariantSubspace
489: @*/
490: PetscErrorCode EPSGetInvariantSubspaceLeft(EPS eps, Vec *v)
491: {
493:   PetscInt       i;

499:   if (!eps->leftvecs) {
500:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
501:   }
502:   if (!eps->W) {
503:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
504:   }
505:   if (!eps->ishermitian && eps->evecsavailable) {
506:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSGetInvariantSubspaceLeft must be called before EPSGetEigenpairLeft,EPSComputeRelativeErrorLeft or EPSComputeResidualNormLeft");
507:   }
508:   for (i=0;i<eps->nconv;i++) {
509:     VecCopy(eps->W[i],v[i]);
510:   }
511:   return(0);
512: }

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

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

522:    Input Parameters:
523: +  eps - eigensolver context 
524: -  i   - index of the solution

526:    Output Parameters:
527: +  eigr - real part of eigenvalue
528: .  eigi - imaginary part of eigenvalue
529: .  Vr   - real part of eigenvector
530: -  Vi   - imaginary part of eigenvector

532:    Notes:
533:    If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is 
534:    configured with complex scalars the eigenvalue is stored 
535:    directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is 
536:    set to zero).

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

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

546:    Level: beginner

548: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSGetEigenvectorLeft(), EPSSolve(), 
549:           EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
550: @*/
551: PetscErrorCode EPSGetEigenpair(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
552: {

557:   if (!eps->eigr || !eps->eigi || !eps->V) {
558:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
559:   }
560:   if (i<0 || i>=eps->nconv) {
561:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
562:   }
563:   EPSGetEigenvalue(eps,i,eigr,eigi);
564:   EPSGetEigenvector(eps,i,Vr,Vi);
565: 
566:   return(0);
567: }

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

574:    Not Collective

576:    Input Parameters:
577: +  eps - eigensolver context 
578: -  i   - index of the solution

580:    Output Parameters:
581: +  eigr - real part of eigenvalue
582: -  eigi - imaginary part of eigenvalue

584:    Notes:
585:    If the eigenvalue is real, then eigi is set to zero. If PETSc is 
586:    configured with complex scalars the eigenvalue is stored 
587:    directly in eigr (eigi is set to zero).

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

593:    Level: beginner

595: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), 
596:           EPSGetEigenpair()
597: @*/
598: PetscErrorCode EPSGetEigenvalue(EPS eps, PetscInt i, PetscScalar *eigr, PetscScalar *eigi)
599: {
600:   PetscInt       k;

604:   if (!eps->eigr || !eps->eigi) {
605:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
606:   }
607:   if (i<0 || i>=eps->nconv) {
608:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
609:   }

611:   if (!eps->perm) k = i;
612:   else k = eps->perm[i];
613: #ifdef PETSC_USE_COMPLEX
614:   if (eigr) *eigr = eps->eigr[k];
615:   if (eigi) *eigi = 0;
616: #else
617:   if (eigr) *eigr = eps->eigr[k];
618:   if (eigi) *eigi = eps->eigi[k];
619: #endif
620: 
621:   return(0);
622: }

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

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

631:    Input Parameters:
632: +  eps - eigensolver context 
633: -  i   - index of the solution

635:    Output Parameters:
636: +  Vr   - real part of eigenvector
637: -  Vi   - imaginary part of eigenvector

639:    Notes:
640:    If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is 
641:    configured with complex scalars the eigenvector is stored 
642:    directly in Vr (Vi is set to zero).

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

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

652:    Level: beginner

654: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), 
655:           EPSGetEigenpair(), EPSGetEigenvectorLeft()
656: @*/
657: PetscErrorCode EPSGetEigenvector(EPS eps, PetscInt i, Vec Vr, Vec Vi)
658: {
660:   PetscInt       k;

664:   if (!eps->V) {
665:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
666:   }
667:   if (i<0 || i>=eps->nconv) {
668:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
669:   }
670:   if (!eps->evecsavailable && (Vr || Vi) ) {
671:     (*eps->ops->computevectors)(eps);
672:   }

674:   if (!eps->perm) k = i;
675:   else k = eps->perm[i];
676: #ifdef PETSC_USE_COMPLEX
677:   if (Vr) { VecCopy(eps->V[k], Vr);  }
678:   if (Vi) { VecSet(Vi,0.0);  }
679: #else
680:   if (eps->eigi[k] > 0) { /* first value of conjugate pair */
681:     if (Vr) { VecCopy(eps->V[k], Vr);  }
682:     if (Vi) { VecCopy(eps->V[k+1], Vi);  }
683:   } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
684:     if (Vr) { VecCopy(eps->V[k-1], Vr);  }
685:     if (Vi) {
686:       VecCopy(eps->V[k], Vi);
687:       VecScale(Vi,-1.0);
688:     }
689:   } else { /* real eigenvalue */
690:     if (Vr) { VecCopy(eps->V[k], Vr);  }
691:     if (Vi) { VecSet(Vi,0.0);  }
692:   }
693: #endif
694: 
695:   return(0);
696: }

700: /*@
701:    EPSGetEigenvectorLeft - Gets the i-th left eigenvector as computed by EPSSolve() 
702:    (only available in two-sided eigensolvers). 

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

706:    Input Parameters:
707: +  eps - eigensolver context 
708: -  i   - index of the solution

710:    Output Parameters:
711: +  Wr   - real part of eigenvector
712: -  Wi   - imaginary part of eigenvector

714:    Notes:
715:    If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is 
716:    configured with complex scalars the eigenvector is stored 
717:    directly in Wr (Wi is set to zero).

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

723:    Level: beginner

725: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), 
726:           EPSGetEigenpair(), EPSGetEigenvector()
727: @*/
728: PetscErrorCode EPSGetEigenvectorLeft(EPS eps, PetscInt i, Vec Wr, Vec Wi)
729: {
731:   PetscInt       k;

735:   if (!eps->leftvecs) {
736:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
737:   }
738:   if (!eps->W) {
739:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
740:   }
741:   if (i<0 || i>=eps->nconv) {
742:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
743:   }
744:   if (!eps->evecsavailable && (Wr || Wi) ) {
745:     (*eps->ops->computevectors)(eps);
746:   }

748:   if (!eps->perm) k = i;
749:   else k = eps->perm[i];
750: #ifdef PETSC_USE_COMPLEX
751:   if (Wr) { VecCopy(eps->W[k], Wr);  }
752:   if (Wi) { VecSet(Wi,0.0);  }
753: #else
754:   if (eps->eigi[k] > 0) { /* first value of conjugate pair */
755:     if (Wr) { VecCopy(eps->W[k], Wr);  }
756:     if (Wi) { VecCopy(eps->W[k+1], Wi);  }
757:   } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
758:     if (Wr) { VecCopy(eps->W[k-1], Wr);  }
759:     if (Wi) {
760:       VecCopy(eps->W[k], Wi);
761:       VecScale(Wi,-1.0);
762:     }
763:   } else { /* real eigenvalue */
764:     if (Wr) { VecCopy(eps->W[k], Wr);  }
765:     if (Wi) { VecSet(Wi,0.0);  }
766:   }
767: #endif
768: 
769:   return(0);
770: }

774: /*@
775:    EPSGetErrorEstimate - Returns the error estimate associated to the i-th 
776:    computed eigenpair.

778:    Not Collective

780:    Input Parameter:
781: +  eps - eigensolver context 
782: -  i   - index of eigenpair

784:    Output Parameter:
785: .  errest - the error estimate

787:    Notes:
788:    This is the error estimate used internally by the eigensolver. The actual
789:    error bound can be computed with EPSComputeRelativeError(). See also the users
790:    manual for details.

792:    Level: advanced

794: .seealso: EPSComputeRelativeError()
795: @*/
796: PetscErrorCode EPSGetErrorEstimate(EPS eps, PetscInt i, PetscReal *errest)
797: {
800:   if (!eps->eigr || !eps->eigi) {
801:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
802:   }
803:   if (i<0 || i>=eps->nconv) {
804:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
805:   }
806:   if (eps->perm) i = eps->perm[i];
807:   if (errest) *errest = eps->errest[i];
808:   return(0);
809: }

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

817:    Not Collective

819:    Input Parameter:
820: +  eps - eigensolver context 
821: -  i   - index of eigenpair

823:    Output Parameter:
824: .  errest - the left error estimate

826:    Notes:
827:    This is the error estimate used internally by the eigensolver. The actual
828:    error bound can be computed with EPSComputeRelativeErrorLeft(). See also the users
829:    manual for details.

831:    Level: advanced

833: .seealso: EPSComputeRelativeErrorLeft()
834: @*/
835: PetscErrorCode EPSGetErrorEstimateLeft(EPS eps, PetscInt i, PetscReal *errest)
836: {
839:   if (!eps->eigr || !eps->eigi) {
840:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
841:   }
842:   if (!eps->leftvecs) {
843:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
844:   }
845:   if (i<0 || i>=eps->nconv) {
846:     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
847:   }
848:   if (eps->perm) i = eps->perm[i];
849:   if (errest) *errest = eps->errest_left[i];
850:   return(0);
851: }

855: /*
856:    EPSComputeResidualNorm_Private - Computes the norm of the residual vector 
857:    associated with an eigenpair.
858: */
859: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps, PetscScalar kr, PetscScalar ki, Vec xr, Vec xi, PetscReal *norm)
860: {
862:   Vec            u, w;
863:   Mat            A, B;
864: #ifndef PETSC_USE_COMPLEX
865:   Vec            v;
866:   PetscReal      ni, nr;
867: #endif
868: 
870:   STGetOperators(eps->OP,&A,&B);
871:   VecDuplicate(eps->V[0],&u);
872:   VecDuplicate(eps->V[0],&w);
873: 
874: #ifndef PETSC_USE_COMPLEX
875:   if (ki == 0 ||
876:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
877: #endif
878:     MatMult(A,xr,u);                             /* u=A*x */
879:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
880:       if (eps->isgeneralized) { MatMult(B,xr,w); }
881:       else { VecCopy(xr,w); }                    /* w=B*x */
882:       VecAXPY(u,-kr,w);                          /* u=A*x-k*B*x */
883:     }
884:     VecNorm(u,NORM_2,norm);
885: #ifndef PETSC_USE_COMPLEX
886:   } else {
887:     VecDuplicate(eps->V[0],&v);
888:     MatMult(A,xr,u);                             /* u=A*xr */
889:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
890:       if (eps->isgeneralized) { MatMult(B,xr,v); }
891:       else { VecCopy(xr,v); }                    /* v=B*xr */
892:       VecAXPY(u,-kr,v);                          /* u=A*xr-kr*B*xr */
893:       if (eps->isgeneralized) { MatMult(B,xi,w); }
894:       else { VecCopy(xi,w); }                    /* w=B*xi */
895:       VecAXPY(u,ki,w);                           /* u=A*xr-kr*B*xr+ki*B*xi */
896:     }
897:     VecNorm(u,NORM_2,&nr);
898:     MatMult(A,xi,u);                             /* u=A*xi */
899:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
900:       VecAXPY(u,-kr,w);                          /* u=A*xi-kr*B*xi */
901:       VecAXPY(u,-ki,v);                          /* u=A*xi-kr*B*xi-ki*B*xr */
902:     }
903:     VecNorm(u,NORM_2,&ni);
904:     *norm = SlepcAbsEigenvalue(nr,ni);
905:     VecDestroy(v);
906:   }
907: #endif

909:   VecDestroy(w);
910:   VecDestroy(u);
911:   return(0);
912: }

916: /*@
917:    EPSComputeResidualNorm - Computes the norm of the residual vector associated with 
918:    the i-th computed eigenpair.

920:    Collective on EPS

922:    Input Parameter:
923: .  eps - the eigensolver context
924: .  i   - the solution index

926:    Output Parameter:
927: .  norm - the residual norm, computed as ||Ax-kBx||_2 where k is the 
928:    eigenvalue and x is the eigenvector. 
929:    If k=0 then the residual norm is computed as ||Ax||_2.

931:    Notes:
932:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
933:    Eigenpairs are indexed according to the ordering criterion established 
934:    with EPSSetWhichEigenpairs().

936:    Level: beginner

938: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
939: @*/
940: PetscErrorCode EPSComputeResidualNorm(EPS eps, PetscInt i, PetscReal *norm)
941: {
943:   Vec            xr, xi;
944:   PetscScalar    kr, ki;
945: 
949:   VecDuplicate(eps->V[0],&xr);
950:   VecDuplicate(eps->V[0],&xi);
951:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
952:   EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,norm);
953:   VecDestroy(xr);
954:   VecDestroy(xi);
955:   return(0);
956: }

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

964:    Collective on EPS

966:    Input Parameter:
967: .  eps - the eigensolver context
968: .  i   - the solution index

970:    Output Parameter:
971: .  norm - the residual norm, computed as ||y'A-ky'B||_2 where k is the 
972:    eigenvalue and y is the left eigenvector. 
973:    If k=0 then the residual norm is computed as ||y'A||_2.

975:    Notes:
976:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
977:    Eigenpairs are indexed according to the ordering criterion established 
978:    with EPSSetWhichEigenpairs().

980:    Level: beginner

982: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
983: @*/
984: PetscErrorCode EPSComputeResidualNormLeft(EPS eps, PetscInt i, PetscReal *norm)
985: {
987:   Vec            u, v, w, xr, xi;
988:   Mat            A, B;
989:   PetscScalar    kr, ki;
990: #ifndef PETSC_USE_COMPLEX
991:   PetscReal      ni, nr;
992: #endif
993: 
996:   if (!eps->leftvecs) {
997:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must request left vectors with EPSSetLeftVectorsWanted");
998:   }
999:   STGetOperators(eps->OP,&A,&B);
1000:   VecDuplicate(eps->W[0],&u);
1001:   VecDuplicate(eps->W[0],&v);
1002:   VecDuplicate(eps->W[0],&w);
1003:   VecDuplicate(eps->W[0],&xr);
1004:   VecDuplicate(eps->W[0],&xi);
1005:   EPSGetEigenvalue(eps,i,&kr,&ki);
1006:   EPSGetEigenvectorLeft(eps,i,xr,xi);

1008: #ifndef PETSC_USE_COMPLEX
1009:   if (ki == 0 ||
1010:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1011: #endif
1012:     MatMultTranspose( A, xr, u );  /* u=A'*x */
1013:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1014:       if (eps->isgeneralized) { MatMultTranspose( B, xr, w );  }
1015:       else { VecCopy( xr, w );  } /* w=B'*x */
1016:       VecAXPY( u, -kr, w);  /* u=A'*x-k*B'*x */
1017:     }
1018:     VecNorm( u, NORM_2, norm);
1019: #ifndef PETSC_USE_COMPLEX
1020:   } else {
1021:     MatMultTranspose( A, xr, u );  /* u=A'*xr */
1022:     if (eps->isgeneralized) { MatMultTranspose( B, xr, v );  }
1023:     else { VecCopy( xr, v );  } /* v=B'*xr */
1024:     VecAXPY( u, -kr, v );  /* u=A'*xr-kr*B'*xr */
1025:     if (eps->isgeneralized) { MatMultTranspose( B, xi, w );  }
1026:     else { VecCopy( xi, w );  } /* w=B'*xi */
1027:     VecAXPY( u, ki, w );  /* u=A'*xr-kr*B'*xr+ki*B'*xi */
1028:     VecNorm( u, NORM_2, &nr );
1029:     MatMultTranspose( A, xi, u );  /* u=A'*xi */
1030:     VecAXPY( u, -kr, w );  /* u=A'*xi-kr*B'*xi */
1031:     VecAXPY( u, -ki, v );  /* u=A'*xi-kr*B'*xi-ki*B'*xr */
1032:     VecNorm( u, NORM_2, &ni );
1033:     *norm = SlepcAbsEigenvalue( nr, ni );
1034:   }
1035: #endif

1037:   VecDestroy(w);
1038:   VecDestroy(v);
1039:   VecDestroy(u);
1040:   VecDestroy(xr);
1041:   VecDestroy(xi);
1042:   return(0);
1043: }

1047: /*
1048:    EPSComputeRelativeError_Private - Computes the relative error bound 
1049:    associated with an eigenpair.
1050: */
1051: PetscErrorCode EPSComputeRelativeError_Private(EPS eps, PetscScalar kr, PetscScalar ki, Vec xr, Vec xi, PetscReal *error)
1052: {
1054:   PetscReal      norm, er;
1055: #ifndef PETSC_USE_COMPLEX
1056:   PetscReal      ei;
1057: #endif
1058: 
1060:   EPSComputeResidualNorm_Private(eps,kr,ki,xr,xi,&norm);

1062: #ifndef PETSC_USE_COMPLEX
1063:   if (ki == 0) {
1064: #endif
1065:     VecNorm(xr,NORM_2,&er);
1066: #ifndef PETSC_USE_COMPLEX
1067:   } else {
1068:     VecNorm(xr,NORM_2,&er);
1069:     VecNorm(xi,NORM_2,&ei);
1070:     er = SlepcAbsEigenvalue(er,ei);
1071:   }
1072: #endif    
1073:   (*eps->conv_func)(eps,kr,ki,norm/er,error,eps->conv_ctx);CHKERRQ(ierr)
1074:   return(0);
1075: }

1079: /*@
1080:    EPSComputeRelativeError - Computes the relative error bound associated 
1081:    with the i-th computed eigenpair.

1083:    Collective on EPS

1085:    Input Parameter:
1086: .  eps - the eigensolver context
1087: .  i   - the solution index

1089:    Output Parameter:
1090: .  error - the relative error bound, computed as ||Ax-kBx||_2/||kx||_2 where 
1091:    k is the eigenvalue and x is the eigenvector. 
1092:    If k=0 the relative error is computed as ||Ax||_2/||x||_2.

1094:    Level: beginner

1096: .seealso: EPSSolve(), EPSComputeResidualNorm(), EPSGetErrorEstimate()
1097: @*/
1098: PetscErrorCode EPSComputeRelativeError(EPS eps, PetscInt i, PetscReal *error)
1099: {
1101:   Vec            xr, xi;
1102:   PetscScalar    kr, ki;
1103: 
1107:   VecDuplicate(eps->V[0],&xr);
1108:   VecDuplicate(eps->V[0],&xi);
1109:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
1110:   EPSComputeRelativeError_Private(eps,kr,ki,xr,xi,error);
1111:   VecDestroy(xr);
1112:   VecDestroy(xi);
1113:   return(0);
1114: }

1118: /*@
1119:    EPSComputeRelativeErrorLeft - Computes the relative error bound associated 
1120:    with the i-th computed eigenvalue and left eigenvector (only available in
1121:    two-sided eigensolvers).

1123:    Collective on EPS

1125:    Input Parameter:
1126: .  eps - the eigensolver context
1127: .  i   - the solution index

1129:    Output Parameter:
1130: .  error - the relative error bound, computed as ||y'A-ky'B||_2/||ky||_2 where 
1131:    k is the eigenvalue and y is the left eigenvector. 
1132:    If k=0 the relative error is computed as ||y'A||_2/||y||_2.

1134:    Level: beginner

1136: .seealso: EPSSolve(), EPSComputeResidualNormLeft(), EPSGetErrorEstimateLeft()
1137: @*/
1138: PetscErrorCode EPSComputeRelativeErrorLeft(EPS eps, PetscInt i, PetscReal *error)
1139: {
1141:   Vec            xr, xi;
1142:   PetscScalar    kr, ki;
1143:   PetscReal      norm, er;
1144: #ifndef PETSC_USE_COMPLEX
1145:   Vec            u;
1146:   PetscReal      ei;
1147: #endif
1148: 
1151:   EPSComputeResidualNormLeft(eps,i,&norm);
1152:   VecDuplicate(eps->W[0],&xr);
1153:   VecDuplicate(eps->W[0],&xi);
1154:   EPSGetEigenvalue(eps,i,&kr,&ki);
1155:   EPSGetEigenvectorLeft(eps,i,xr,xi);

1157: #ifndef PETSC_USE_COMPLEX
1158:   if (ki == 0 ||
1159:     PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1160: #endif
1161:     VecNorm(xr, NORM_2, &er);
1162:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1163:       *error =  norm / (PetscAbsScalar(kr) * er);
1164:     } else {
1165:       *error = norm / er;
1166:     }
1167: #ifndef PETSC_USE_COMPLEX
1168:   } else {
1169:     VecDuplicate(xi, &u);
1170:     VecCopy(xi, u);
1171:     VecAXPBY(u, kr, -ki, xr);
1172:     VecNorm(u, NORM_2, &er);
1173:     VecAXPBY(xi, kr, ki, xr);
1174:     VecNorm(xi, NORM_2, &ei);
1175:     VecDestroy(u);
1176:     *error = norm / SlepcAbsEigenvalue(er, ei);
1177:   }
1178: #endif    
1179: 
1180:   VecDestroy(xr);
1181:   VecDestroy(xi);
1182:   return(0);
1183: }

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

1189: /*@
1190:    EPSSortEigenvalues - Sorts a list of eigenvalues according to the criterion 
1191:    specified via EPSSetWhichEigenpairs().

1193:    Not Collective

1195:    Input Parameters:
1196: +  eps   - the eigensolver context
1197: .  n     - number of eigenvalues in the list
1198: .  eigr  - pointer to the array containing the eigenvalues
1199: -  eigi  - imaginary part of the eigenvalues (only when using real numbers)

1201:    Output Parameter:
1202: .  perm  - resulting permutation

1204:    Note:
1205:    The result is a list of indices in the original eigenvalue array 
1206:    corresponding to the first nev eigenvalues sorted in the specified
1207:    criterion.

1209:    Level: developer

1211: .seealso: EPSSortEigenvaluesReal(), EPSSetWhichEigenpairs()
1212: @*/
1213: PetscErrorCode EPSSortEigenvalues(EPS eps,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,PetscInt *perm)
1214: {
1216:   PetscScalar    re,im;
1217:   PetscInt       i,j,result,tmp;

1220:   for (i=0; i<n; i++) { perm[i] = i; }
1221:   /* insertion sort */
1222:   for (i=n-1; i>=0; i--) {
1223:     re = eigr[perm[i]];
1224:     im = eigi[perm[i]];
1225:     j = i + 1;
1226: #ifndef PETSC_USE_COMPLEX
1227:     if (im != 0) {
1228:       /* complex eigenvalue */
1229:       i--;
1230:       im = eigi[perm[i]];
1231:     }
1232: #endif
1233:     while (j<n) {
1234:       EPSCompareEigenvalues(eps,re,im,eigr[perm[j]],eigi[perm[j]],&result);
1235:       if (result < 0) break;
1236: #ifndef PETSC_USE_COMPLEX
1237:       /* keep together every complex conjugated eigenpair */
1238:       if (im == 0) {
1239:         if (eigi[perm[j]] == 0) {
1240: #endif
1241:           tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = tmp;
1242:           j++;
1243: #ifndef PETSC_USE_COMPLEX
1244:         } else {
1245:           tmp = perm[j-1]; perm[j-1] = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp;
1246:           j+=2;
1247:         }
1248:       } else {
1249:         if (eigi[perm[j]] == 0) {
1250:           tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = perm[j-1]; perm[j-1] = tmp;
1251:           j++;
1252:         } else {
1253:           tmp = perm[j-2]; perm[j-2] = perm[j]; perm[j] = tmp;
1254:           tmp = perm[j-1]; perm[j-1] = perm[j+1]; perm[j+1] = tmp;
1255:           j+=2;
1256:         }
1257:       }
1258: #endif
1259:     }
1260:   }
1261:   return(0);
1262: }

1266: /*@
1267:    EPSSortEigenvaluesReal - Sorts a list of eigenvalues according to a certain
1268:    criterion (version for real eigenvalues only).

1270:    Not Collective

1272:    Input Parameters:
1273: +  eps   - the eigensolver context
1274: .  n     - number of eigenvalue in the list
1275: -  eig   - pointer to the array containing the eigenvalues (real)

1277:    Output Parameter:
1278: .  perm  - resulting permutation

1280:    Note:
1281:    The result is a list of indices in the original eigenvalue array 
1282:    corresponding to the first nev eigenvalues sorted in the specified
1283:    criterion.

1285:    Level: developer

1287: .seealso: EPSSortEigenvalues(), EPSSetWhichEigenpairs(), EPSCompareEigenvalues()
1288: @*/
1289: PetscErrorCode EPSSortEigenvaluesReal(EPS eps,PetscInt n,PetscReal *eig,PetscInt *perm)
1290: {
1292:   PetscScalar    re;
1293:   PetscInt       i,j,result,tmp;

1296:   for (i=0; i<n; i++) { perm[i] = i; }
1297:   /* insertion sort */
1298:   for (i=1; i<n; i++) {
1299:     re = eig[perm[i]];
1300:     j = i-1;
1301:     EPSCompareEigenvalues(eps,re,0.0,eig[perm[j]],0.0,&result);
1302:     while (result<=0 && j>=0) {
1303:       tmp = perm[j]; perm[j] = perm[j+1]; perm[j+1] = tmp; j--;
1304:       if (j>=0) {
1305:         EPSCompareEigenvalues(eps,re,0.0,eig[perm[j]],0.0,&result);
1306:       }
1307:     }
1308:   }
1309:   return(0);
1310: }

1314: /*@
1315:    EPSCompareEigenvalues - Compares two (possibly complex) eigenvalues according
1316:    to a certain criterion.

1318:    Not Collective

1320:    Input Parameters:
1321: +  eps   - the eigensolver context
1322: .  ar     - real part of the 1st eigenvalue
1323: .  ai     - imaginary part of the 1st eigenvalue
1324: .  br     - real part of the 2nd eigenvalue
1325: -  bi     - imaginary part of the 2nd eigenvalue

1327:    Output Parameter:
1328: .  res    - result of comparison

1330:    Notes:
1331:    The returning parameter 'res' can be:
1332: +  negative - if the 1st eigenvalue is preferred to the 2st one
1333: .  zero     - if both eigenvalues are equally preferred
1334: -  positive - if the 2st eigenvalue is preferred to the 1st one

1336:    The criterion of comparison is related to the 'which' parameter set with
1337:    EPSSetWhichEigenpairs().

1339:    Level: developer

1341: .seealso: EPSSortEigenvalues(), EPSSetWhichEigenpairs()
1342: @*/
1343: PetscErrorCode EPSCompareEigenvalues(EPS eps,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *result)
1344: {
1346:   PetscReal      a,b;

1349:   switch(eps->which) {
1350:     case EPS_WHICH_USER:
1351:       if (!eps->which_func) SETERRQ(1,"Undefined eigenvalue comparison function");
1352:       (*eps->which_func)(eps,ar,ai,br,bi,result,eps->which_ctx);
1353:       break;
1354:     case EPS_LARGEST_MAGNITUDE:
1355:     case EPS_SMALLEST_MAGNITUDE:
1356:       a = SlepcAbsEigenvalue(ar,ai);
1357:       b = SlepcAbsEigenvalue(br,bi);
1358:       break;
1359:     case EPS_LARGEST_REAL:
1360:     case EPS_SMALLEST_REAL:
1361:       a = PetscRealPart(ar);
1362:       b = PetscRealPart(br);
1363:       break;
1364:     case EPS_LARGEST_IMAGINARY:
1365:     case EPS_SMALLEST_IMAGINARY:
1366: #if defined(PETSC_USE_COMPLEX)
1367:       a = PetscImaginaryPart(ar);
1368:       b = PetscImaginaryPart(br);
1369: #else
1370:       a = PetscAbsReal(ai);
1371:       b = PetscAbsReal(bi);
1372: #endif
1373:       break;
1374:     case EPS_TARGET_MAGNITUDE:
1375:       /* complex target only allowed if scalartype=complex */
1376:       a = SlepcAbsEigenvalue(ar-eps->target,ai);
1377:       b = SlepcAbsEigenvalue(br-eps->target,bi);
1378:       break;
1379:     case EPS_TARGET_REAL:
1380:       a = PetscAbsReal(PetscRealPart(ar-eps->target));
1381:       b = PetscAbsReal(PetscRealPart(br-eps->target));
1382:       break;
1383:     case EPS_TARGET_IMAGINARY:
1384: #if !defined(PETSC_USE_COMPLEX)
1385:       /* complex target only allowed if scalartype=complex */
1386:       a = PetscAbsReal(ai);
1387:       b = PetscAbsReal(bi);
1388: #else
1389:       a = PetscAbsReal(PetscImaginaryPart(ar-eps->target));
1390:       b = PetscAbsReal(PetscImaginaryPart(br-eps->target));
1391: #endif
1392:       break;
1393:     default: SETERRQ(1,"Wrong value of which");
1394:   }
1395:   switch(eps->which) {
1396:     case EPS_WHICH_USER:
1397:       break;
1398:     case EPS_LARGEST_MAGNITUDE:
1399:     case EPS_LARGEST_REAL:
1400:     case EPS_LARGEST_IMAGINARY:
1401:       if (a<b) *result = 1;
1402:       else if (a>b) *result = -1;
1403:       else *result = 0;
1404:       break;
1405:     default:
1406:       if (a>b) *result = 1;
1407:       else if (a<b) *result = -1;
1408:       else *result = 0;
1409:   }
1410:   return(0);
1411: }

1415: /*@
1416:    EPSGetStartVector - Gets a suitable vector to be used as the starting vector
1417:    for the recurrence that builds the right subspace.

1419:    Collective on EPS and Vec

1421:    Input Parameters:
1422: +  eps - the eigensolver context
1423: -  i   - iteration number

1425:    Output Parameters:
1426: +  vec - the start vector
1427: -  breakdown - flag indicating that a breakdown has occurred

1429:    Notes:
1430:    The start vector is computed from another vector: for the first step (i=0),
1431:    the first initial vector is used (see EPSSetInitialSpace()); otherwise a random
1432:    vector is created. Then this vector is forced to be in the range of OP (only
1433:    for generalized definite problems) and orthonormalized with respect to all
1434:    V-vectors up to i-1.

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

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

1443:    Level: developer

1445: .seealso: EPSSetInitialSpace()
1446: @*/
1447: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,Vec vec,PetscTruth *breakdown)
1448: {
1450:   PetscReal      norm;
1451:   PetscTruth     lindep;
1452:   Vec            w;
1453: 

1458:   VecDuplicate(eps->V[0],&w);

1460:   /* For the first step, use the first initial vector, otherwise a random one */
1461:   if (i==0 && eps->nini>0) {
1462:     VecCopy(eps->V[0],w);
1463:   } else {
1464:     SlepcVecSetRandom(w,eps->rand);
1465:   }

1467:   /* Force the vector to be in the range of OP for definite generalized problems */
1468:   if (eps->ispositive) {
1469:     STApply(eps->OP,w,vec);
1470:   } else {
1471:     VecCopy(w,vec);
1472:   }

1474:   /* Orthonormalize the vector with respect to previous vectors */
1475:   IPOrthogonalize(eps->ip,eps->nds,eps->DS,i,PETSC_NULL,eps->V,vec,PETSC_NULL,&norm,&lindep);
1476:   if (breakdown) *breakdown = lindep;
1477:   else if (lindep || norm == 0.0) {
1478:     if (i==0) { SETERRQ(1,"Initial vector is zero or belongs to the deflation space"); }
1479:     else { SETERRQ(1,"Unable to generate more start vectors"); }
1480:   }
1481:   VecScale(vec,1.0/norm);

1483:   VecDestroy(w);
1484:   return(0);
1485: }

1489: /*@
1490:    EPSGetStartVectorLeft - Gets a suitable vector to be used as the starting vector
1491:    in the recurrence that builds the left subspace (in methods that work with two
1492:    subspaces).

1494:    Collective on EPS and Vec

1496:    Input Parameters:
1497: +  eps - the eigensolver context
1498: -  i   - iteration number

1500:    Output Parameter:
1501: +  vec - the start vector
1502: -  breakdown - flag indicating that a breakdown has occurred

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

1510:    The flag breakdown is set to true if i>0 and the vector is linearly dependent
1511:    with respect to the W-vectors.

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

1516:    Level: developer

1518: .seealso: EPSSetInitialSpaceLeft()

1520: @*/
1521: PetscErrorCode EPSGetStartVectorLeft(EPS eps,PetscInt i,Vec vec,PetscTruth *breakdown)
1522: {
1524:   PetscReal      norm;
1525:   PetscTruth     lindep;
1526:   Vec            w;
1527: 

1532:   VecDuplicate(eps->W[0],&w);

1534:   /* For the first step, use the first initial left vector, otherwise a random one */
1535:   if (i==0 && eps->ninil>0) {
1536:     VecCopy(eps->W[0],w);
1537:   } else {
1538:     SlepcVecSetRandom(w,eps->rand);
1539:   }

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

1544:   /* Orthonormalize the vector with respect to previous vectors */
1545:   IPOrthogonalize(eps->ip,0,PETSC_NULL,i,PETSC_NULL,eps->W,vec,PETSC_NULL,&norm,&lindep);
1546:   if (breakdown) *breakdown = lindep;
1547:   else if (lindep || norm == 0.0) {
1548:     if (i==0) { SETERRQ(1,"Left initial vector is zero"); }
1549:     else { SETERRQ(1,"Unable to generate more left start vectors"); }
1550:   }
1551:   VecScale(vec,1/norm);

1553:   VecDestroy(w);
1554:   return(0);
1555: }