Actual source code: epssolve.c

slepc-3.17.0 2022-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    EPS routines related to the solution process
 12: */

 14: #include <slepc/private/epsimpl.h>
 15: #include <slepc/private/bvimpl.h>
 16: #include <petscdraw.h>

 18: PetscErrorCode EPSComputeVectors(EPS eps)
 19: {
 20:   EPSCheckSolved(eps,1);
 21:   if (eps->state==EPS_STATE_SOLVED && eps->ops->computevectors) (*eps->ops->computevectors)(eps);
 22:   eps->state = EPS_STATE_EIGENVECTORS;
 23:   PetscFunctionReturn(0);
 24: }

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

 28: static PetscErrorCode EPSComputeValues(EPS eps)
 29: {
 30:   PetscBool      injective,iscomp,isfilter;
 31:   PetscInt       i,n,aux,nconv0;
 32:   Mat            A,B=NULL,G,Z;

 34:   switch (eps->categ) {
 35:     case EPS_CATEGORY_KRYLOV:
 36:     case EPS_CATEGORY_OTHER:
 37:       STIsInjective(eps->st,&injective);
 38:       if (injective) {
 39:         /* one-to-one mapping: backtransform eigenvalues */
 41:         (*eps->ops->backtransform)(eps);
 42:       } else {
 43:         /* compute eigenvalues from Rayleigh quotient */
 44:         DSGetDimensions(eps->ds,&n,NULL,NULL,NULL);
 45:         if (!n) break;
 46:         EPSGetOperators(eps,&A,&B);
 47:         BVSetActiveColumns(eps->V,0,n);
 48:         DSGetCompact(eps->ds,&iscomp);
 49:         DSSetCompact(eps->ds,PETSC_FALSE);
 50:         DSGetMat(eps->ds,DS_MAT_A,&G);
 51:         BVMatProject(eps->V,A,eps->V,G);
 52:         DSRestoreMat(eps->ds,DS_MAT_A,&G);
 53:         if (B) {
 54:           DSGetMat(eps->ds,DS_MAT_B,&G);
 55:           BVMatProject(eps->V,B,eps->V,G);
 56:           DSRestoreMat(eps->ds,DS_MAT_A,&G);
 57:         }
 58:         DSSolve(eps->ds,eps->eigr,eps->eigi);
 59:         DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
 60:         DSSynchronize(eps->ds,eps->eigr,eps->eigi);
 61:         DSSetCompact(eps->ds,iscomp);
 62:         if (eps->ishermitian && (!eps->isgeneralized || eps->ispositive)) { /* V = V * Z */
 63:           DSVectors(eps->ds,DS_MAT_X,NULL,NULL);
 64:           DSGetMat(eps->ds,DS_MAT_X,&Z);
 65:           BVMultInPlace(eps->V,Z,0,n);
 66:           MatDestroy(&Z);
 67:         }
 68:         /* in case of STFILTER discard computed eigenvalues that lie outside the wanted interval */
 69:         PetscObjectTypeCompare((PetscObject)eps->st,STFILTER,&isfilter);
 70:         if (isfilter) {
 71:           nconv0 = eps->nconv;
 72:           for (i=0;i<eps->nconv;i++) {
 73:             if (PetscRealPart(eps->eigr[eps->perm[i]])<eps->inta || PetscRealPart(eps->eigr[eps->perm[i]])>eps->intb) {
 74:               eps->nconv--;
 75:               if (i<eps->nconv) { SWAP(eps->perm[i],eps->perm[eps->nconv],aux); i--; }
 76:             }
 77:           }
 78:           if (nconv0>eps->nconv) PetscInfo(eps,"Discarded %" PetscInt_FMT " computed eigenvalues lying outside the interval\n",nconv0-eps->nconv);
 79:         }
 80:       }
 81:       break;
 82:     case EPS_CATEGORY_PRECOND:
 83:     case EPS_CATEGORY_CONTOUR:
 84:       /* eigenvalues already available as an output of the solver */
 85:       break;
 86:   }
 87:   PetscFunctionReturn(0);
 88: }

 90: /*@
 91:    EPSSolve - Solves the eigensystem.

 93:    Collective on eps

 95:    Input Parameter:
 96: .  eps - eigensolver context obtained from EPSCreate()

 98:    Options Database Keys:
 99: +  -eps_view - print information about the solver used
100: .  -eps_view_mat0 - view the first matrix (A)
101: .  -eps_view_mat1 - view the second matrix (B)
102: .  -eps_view_vectors - view the computed eigenvectors
103: .  -eps_view_values - view the computed eigenvalues
104: .  -eps_converged_reason - print reason for convergence, and number of iterations
105: .  -eps_error_absolute - print absolute errors of each eigenpair
106: .  -eps_error_relative - print relative errors of each eigenpair
107: -  -eps_error_backward - print backward errors of each eigenpair

109:    Notes:
110:    All the command-line options listed above admit an optional argument specifying
111:    the viewer type and options. For instance, use '-eps_view_mat0 binary:amatrix.bin'
112:    to save the A matrix to a binary file, '-eps_view_values draw' to draw the computed
113:    eigenvalues graphically, or '-eps_error_relative :myerr.m:ascii_matlab' to save
114:    the errors in a file that can be executed in Matlab.

116:    Level: beginner

118: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
119: @*/
120: PetscErrorCode EPSSolve(EPS eps)
121: {
122:   PetscInt       i;
123:   STMatMode      matmode;
124:   Mat            A,B;

127:   if (eps->state>=EPS_STATE_SOLVED) PetscFunctionReturn(0);
128:   PetscLogEventBegin(EPS_Solve,eps,0,0,0);

130:   /* Call setup */
131:   EPSSetUp(eps);
132:   eps->nconv = 0;
133:   eps->its   = 0;
134:   for (i=0;i<eps->ncv;i++) {
135:     eps->eigr[i]   = 0.0;
136:     eps->eigi[i]   = 0.0;
137:     eps->errest[i] = 0.0;
138:     eps->perm[i]   = i;
139:   }
140:   EPSViewFromOptions(eps,NULL,"-eps_view_pre");
141:   RGViewFromOptions(eps->rg,NULL,"-rg_view");

143:   /* Call solver */
144:   (*eps->ops->solve)(eps);
146:   eps->state = EPS_STATE_SOLVED;

148:   /* Only the first nconv columns contain useful information (except in CISS) */
149:   BVSetActiveColumns(eps->V,0,eps->nconv);
150:   if (eps->twosided) BVSetActiveColumns(eps->W,0,eps->nconv);

152:   /* If inplace, purify eigenvectors before reverting operator */
153:   STGetMatMode(eps->st,&matmode);
154:   if (matmode == ST_MATMODE_INPLACE && eps->ispositive) EPSComputeVectors(eps);
155:   STPostSolve(eps->st);

157:   /* Map eigenvalues back to the original problem if appropriate */
158:   EPSComputeValues(eps);

160: #if !defined(PETSC_USE_COMPLEX)
161:   /* Reorder conjugate eigenvalues (positive imaginary first) */
162:   for (i=0;i<eps->nconv-1;i++) {
163:     if (eps->eigi[i] != 0) {
164:       if (eps->eigi[i] < 0) {
165:         eps->eigi[i] = -eps->eigi[i];
166:         eps->eigi[i+1] = -eps->eigi[i+1];
167:         /* the next correction only works with eigenvectors */
168:         EPSComputeVectors(eps);
169:         BVScaleColumn(eps->V,i+1,-1.0);
170:       }
171:       i++;
172:     }
173:   }
174: #endif

176:   /* Sort eigenvalues according to eps->which parameter */
177:   SlepcSortEigenvalues(eps->sc,eps->nconv,eps->eigr,eps->eigi,eps->perm);
178:   PetscLogEventEnd(EPS_Solve,eps,0,0,0);

180:   /* Various viewers */
181:   EPSViewFromOptions(eps,NULL,"-eps_view");
182:   EPSConvergedReasonViewFromOptions(eps);
183:   EPSErrorViewFromOptions(eps);
184:   EPSValuesViewFromOptions(eps);
185:   EPSVectorsViewFromOptions(eps);
186:   EPSGetOperators(eps,&A,&B);
187:   MatViewFromOptions(A,(PetscObject)eps,"-eps_view_mat0");
188:   if (eps->isgeneralized) MatViewFromOptions(B,(PetscObject)eps,"-eps_view_mat1");

190:   /* Remove deflation and initial subspaces */
191:   if (eps->nds) {
192:     BVSetNumConstraints(eps->V,0);
193:     eps->nds = 0;
194:   }
195:   eps->nini = 0;
196:   PetscFunctionReturn(0);
197: }

199: /*@
200:    EPSGetIterationNumber - Gets the current iteration number. If the
201:    call to EPSSolve() is complete, then it returns the number of iterations
202:    carried out by the solution method.

204:    Not Collective

206:    Input Parameter:
207: .  eps - the eigensolver context

209:    Output Parameter:
210: .  its - number of iterations

212:    Note:
213:    During the i-th iteration this call returns i-1. If EPSSolve() is
214:    complete, then parameter "its" contains either the iteration number at
215:    which convergence was successfully reached, or failure was detected.
216:    Call EPSGetConvergedReason() to determine if the solver converged or
217:    failed and why.

219:    Level: intermediate

221: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
222: @*/
223: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
224: {
227:   *its = eps->its;
228:   PetscFunctionReturn(0);
229: }

231: /*@
232:    EPSGetConverged - Gets the number of converged eigenpairs.

234:    Not Collective

236:    Input Parameter:
237: .  eps - the eigensolver context

239:    Output Parameter:
240: .  nconv - number of converged eigenpairs

242:    Note:
243:    This function should be called after EPSSolve() has finished.

245:    Level: beginner

247: .seealso: EPSSetDimensions(), EPSSolve(), EPSGetEigenpair()
248: @*/
249: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
250: {
253:   EPSCheckSolved(eps,1);
254:   *nconv = eps->nconv;
255:   PetscFunctionReturn(0);
256: }

258: /*@
259:    EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
260:    stopped.

262:    Not Collective

264:    Input Parameter:
265: .  eps - the eigensolver context

267:    Output Parameter:
268: .  reason - negative value indicates diverged, positive value converged

270:    Options Database Key:
271: .  -eps_converged_reason - print the reason to a viewer

273:    Notes:
274:    Possible values for reason are
275: +  EPS_CONVERGED_TOL - converged up to tolerance
276: .  EPS_CONVERGED_USER - converged due to a user-defined condition
277: .  EPS_DIVERGED_ITS - required more than max_it iterations to reach convergence
278: .  EPS_DIVERGED_BREAKDOWN - generic breakdown in method
279: -  EPS_DIVERGED_SYMMETRY_LOST - pseudo-Lanczos was not able to keep symmetry

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

283:    Level: intermediate

285: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
286: @*/
287: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
288: {
291:   EPSCheckSolved(eps,1);
292:   *reason = eps->reason;
293:   PetscFunctionReturn(0);
294: }

296: /*@
297:    EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
298:    subspace.

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

302:    Input Parameter:
303: .  eps - the eigensolver context

305:    Output Parameter:
306: .  v - an array of vectors

308:    Notes:
309:    This function should be called after EPSSolve() has finished.

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

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

321:    Level: intermediate

323: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve()
324: @*/
325: PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec v[])
326: {
327:   PetscInt       i;
328:   BV             V=eps->V;
329:   Vec            w;

334:   EPSCheckSolved(eps,1);
336:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
337:     BVDuplicateResize(eps->V,eps->nconv,&V);
338:     BVSetActiveColumns(eps->V,0,eps->nconv);
339:     BVCopy(eps->V,V);
340:     for (i=0;i<eps->nconv;i++) {
341:       BVGetColumn(V,i,&w);
342:       VecPointwiseDivide(w,w,eps->D);
343:       BVRestoreColumn(V,i,&w);
344:     }
345:     BVOrthogonalize(V,NULL);
346:   }
347:   for (i=0;i<eps->nconv;i++) BVCopyVec(V,i,v[i]);
348:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) BVDestroy(&V);
349:   PetscFunctionReturn(0);
350: }

352: /*@C
353:    EPSGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
354:    EPSSolve(). The solution consists in both the eigenvalue and the eigenvector.

356:    Logically Collective on eps

358:    Input Parameters:
359: +  eps - eigensolver context
360: -  i   - index of the solution

362:    Output Parameters:
363: +  eigr - real part of eigenvalue
364: .  eigi - imaginary part of eigenvalue
365: .  Vr   - real part of eigenvector
366: -  Vi   - imaginary part of eigenvector

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

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

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

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

386:    Level: beginner

388: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSGetLeftEigenvector(), EPSSolve(),
389:           EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
390: @*/
391: PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
392: {
395:   EPSCheckSolved(eps,1);
398:   EPSGetEigenvalue(eps,i,eigr,eigi);
399:   if (Vr || Vi) EPSGetEigenvector(eps,i,Vr,Vi);
400:   PetscFunctionReturn(0);
401: }

403: /*@C
404:    EPSGetEigenvalue - Gets the i-th eigenvalue as computed by EPSSolve().

406:    Not Collective

408:    Input Parameters:
409: +  eps - eigensolver context
410: -  i   - index of the solution

412:    Output Parameters:
413: +  eigr - real part of eigenvalue
414: -  eigi - imaginary part of eigenvalue

416:    Notes:
417:    If the eigenvalue is real, then eigi is set to zero. If PETSc is
418:    configured with complex scalars the eigenvalue is stored
419:    directly in eigr (eigi is set to zero).

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

425:    Level: beginner

427: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
428: @*/
429: PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi)
430: {
431:   PetscInt k;

434:   EPSCheckSolved(eps,1);
437:   k = eps->perm[i];
438: #if defined(PETSC_USE_COMPLEX)
439:   if (eigr) *eigr = eps->eigr[k];
440:   if (eigi) *eigi = 0;
441: #else
442:   if (eigr) *eigr = eps->eigr[k];
443:   if (eigi) *eigi = eps->eigi[k];
444: #endif
445:   PetscFunctionReturn(0);
446: }

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

451:    Logically Collective on eps

453:    Input Parameters:
454: +  eps - eigensolver context
455: -  i   - index of the solution

457:    Output Parameters:
458: +  Vr   - real part of eigenvector
459: -  Vi   - imaginary part of eigenvector

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

465:    If the corresponding eigenvalue is real, then Vi is set to zero. If PETSc is
466:    configured with complex scalars the eigenvector is stored
467:    directly in Vr (Vi is set to zero). In any case, the user can pass NULL in Vr
468:    or Vi if one of them is not required.

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

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

478:    Level: beginner

480: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair(), EPSGetLeftEigenvector()
481: @*/
482: PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi)
483: {
484:   PetscInt       k;

490:   EPSCheckSolved(eps,1);
493:   EPSComputeVectors(eps);
494:   k = eps->perm[i];
495:   BV_GetEigenvector(eps->V,k,eps->eigi[k],Vr,Vi);
496:   PetscFunctionReturn(0);
497: }

499: /*@
500:    EPSGetLeftEigenvector - Gets the i-th left eigenvector as computed by EPSSolve().

502:    Logically Collective on eps

504:    Input Parameters:
505: +  eps - eigensolver context
506: -  i   - index of the solution

508:    Output Parameters:
509: +  Wr   - real part of left eigenvector
510: -  Wi   - imaginary part of left eigenvector

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

516:    If the corresponding eigenvalue is real, then Wi is set to zero. If PETSc is
517:    configured with complex scalars the eigenvector is stored directly in Wr
518:    (Wi is set to zero). In any case, the user can pass NULL in Wr or Wi if
519:    one of them is not required.

521:    The index i should be a value between 0 and nconv-1 (see EPSGetConverged()).
522:    Eigensolutions are indexed according to the ordering criterion established
523:    with EPSSetWhichEigenpairs().

525:    Left eigenvectors are available only if the twosided flag was set, see
526:    EPSSetTwoSided().

528:    Level: intermediate

530: .seealso: EPSGetEigenvector(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSSetTwoSided()
531: @*/
532: PetscErrorCode EPSGetLeftEigenvector(EPS eps,PetscInt i,Vec Wr,Vec Wi)
533: {
534:   PetscInt       k;

540:   EPSCheckSolved(eps,1);
544:   EPSComputeVectors(eps);
545:   k = eps->perm[i];
546:   BV_GetEigenvector(eps->W,k,eps->eigi[k],Wr,Wi);
547:   PetscFunctionReturn(0);
548: }

550: /*@
551:    EPSGetErrorEstimate - Returns the error estimate associated to the i-th
552:    computed eigenpair.

554:    Not Collective

556:    Input Parameters:
557: +  eps - eigensolver context
558: -  i   - index of eigenpair

560:    Output Parameter:
561: .  errest - the error estimate

563:    Notes:
564:    This is the error estimate used internally by the eigensolver. The actual
565:    error bound can be computed with EPSComputeError(). See also the users
566:    manual for details.

568:    Level: advanced

570: .seealso: EPSComputeError()
571: @*/
572: PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)
573: {
576:   EPSCheckSolved(eps,1);
579:   *errest = eps->errest[eps->perm[i]];
580:   PetscFunctionReturn(0);
581: }

583: /*
584:    EPSComputeResidualNorm_Private - Computes the norm of the residual vector
585:    associated with an eigenpair.

587:    Input Parameters:
588:      trans - whether A' must be used instead of A
589:      kr,ki - eigenvalue
590:      xr,xi - eigenvector
591:      z     - three work vectors (the second one not referenced in complex scalars)
592: */
593: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
594: {
595:   PetscInt       nmat;
596:   Mat            A,B;
597:   Vec            u,w;
598:   PetscScalar    alpha;
599: #if !defined(PETSC_USE_COMPLEX)
600:   Vec            v;
601:   PetscReal      ni,nr;
602: #endif
603:   PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultHermitianTranspose: MatMult;

605:   u = z[0]; w = z[2];
606:   STGetNumMatrices(eps->st,&nmat);
607:   STGetMatrix(eps->st,0,&A);
608:   if (nmat>1) STGetMatrix(eps->st,1,&B);

610: #if !defined(PETSC_USE_COMPLEX)
611:   v = z[1];
612:   if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
613: #endif
614:     (*matmult)(A,xr,u);                          /* u=A*x */
615:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
616:       if (nmat>1) (*matmult)(B,xr,w);
617:       else VecCopy(xr,w);                        /* w=B*x */
618:       alpha = trans? -PetscConj(kr): -kr;
619:       VecAXPY(u,alpha,w);                        /* u=A*x-k*B*x */
620:     }
621:     VecNorm(u,NORM_2,norm);
622: #if !defined(PETSC_USE_COMPLEX)
623:   } else {
624:     (*matmult)(A,xr,u);                          /* u=A*xr */
625:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
626:       if (nmat>1) (*matmult)(B,xr,v);
627:       else VecCopy(xr,v);                        /* v=B*xr */
628:       VecAXPY(u,-kr,v);                          /* u=A*xr-kr*B*xr */
629:       if (nmat>1) (*matmult)(B,xi,w);
630:       else VecCopy(xi,w);                        /* w=B*xi */
631:       VecAXPY(u,trans?-ki:ki,w);                 /* u=A*xr-kr*B*xr+ki*B*xi */
632:     }
633:     VecNorm(u,NORM_2,&nr);
634:     (*matmult)(A,xi,u);                          /* u=A*xi */
635:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
636:       VecAXPY(u,-kr,w);                          /* u=A*xi-kr*B*xi */
637:       VecAXPY(u,trans?ki:-ki,v);                 /* u=A*xi-kr*B*xi-ki*B*xr */
638:     }
639:     VecNorm(u,NORM_2,&ni);
640:     *norm = SlepcAbsEigenvalue(nr,ni);
641:   }
642: #endif
643:   PetscFunctionReturn(0);
644: }

646: /*@
647:    EPSComputeError - Computes the error (based on the residual norm) associated
648:    with the i-th computed eigenpair.

650:    Collective on eps

652:    Input Parameters:
653: +  eps  - the eigensolver context
654: .  i    - the solution index
655: -  type - the type of error to compute

657:    Output Parameter:
658: .  error - the error

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

664:    Level: beginner

666: .seealso: EPSErrorType, EPSSolve(), EPSGetErrorEstimate()
667: @*/
668: PetscErrorCode EPSComputeError(EPS eps,PetscInt i,EPSErrorType type,PetscReal *error)
669: {
670:   Mat            A,B;
671:   Vec            xr,xi,w[3];
672:   PetscReal      t,vecnorm=1.0,errorl;
673:   PetscScalar    kr,ki;
674:   PetscBool      flg;

680:   EPSCheckSolved(eps,1);

682:   /* allocate work vectors */
683: #if defined(PETSC_USE_COMPLEX)
684:   EPSSetWorkVecs(eps,3);
685:   xi   = NULL;
686:   w[1] = NULL;
687: #else
688:   EPSSetWorkVecs(eps,5);
689:   xi   = eps->work[3];
690:   w[1] = eps->work[4];
691: #endif
692:   xr   = eps->work[0];
693:   w[0] = eps->work[1];
694:   w[2] = eps->work[2];

696:   /* compute residual norm */
697:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
698:   EPSComputeResidualNorm_Private(eps,PETSC_FALSE,kr,ki,xr,xi,w,error);

700:   /* compute 2-norm of eigenvector */
701:   if (eps->problem_type==EPS_GHEP) VecNorm(xr,NORM_2,&vecnorm);

703:   /* if two-sided, compute left residual norm and take the maximum */
704:   if (eps->twosided) {
705:     EPSGetLeftEigenvector(eps,i,xr,xi);
706:     EPSComputeResidualNorm_Private(eps,PETSC_TRUE,kr,ki,xr,xi,w,&errorl);
707:     *error = PetscMax(*error,errorl);
708:   }

710:   /* compute error */
711:   switch (type) {
712:     case EPS_ERROR_ABSOLUTE:
713:       break;
714:     case EPS_ERROR_RELATIVE:
715:       *error /= SlepcAbsEigenvalue(kr,ki)*vecnorm;
716:       break;
717:     case EPS_ERROR_BACKWARD:
718:       /* initialization of matrix norms */
719:       if (!eps->nrma) {
720:         STGetMatrix(eps->st,0,&A);
721:         MatHasOperation(A,MATOP_NORM,&flg);
723:         MatNorm(A,NORM_INFINITY,&eps->nrma);
724:       }
725:       if (eps->isgeneralized) {
726:         if (!eps->nrmb) {
727:           STGetMatrix(eps->st,1,&B);
728:           MatHasOperation(B,MATOP_NORM,&flg);
730:           MatNorm(B,NORM_INFINITY,&eps->nrmb);
731:         }
732:       } else eps->nrmb = 1.0;
733:       t = SlepcAbsEigenvalue(kr,ki);
734:       *error /= (eps->nrma+t*eps->nrmb)*vecnorm;
735:       break;
736:     default:
737:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
738:   }
739:   PetscFunctionReturn(0);
740: }

742: /*
743:    EPSGetStartVector - Generate a suitable vector to be used as the starting vector
744:    for the recurrence that builds the right subspace.

746:    Collective on eps

748:    Input Parameters:
749: +  eps - the eigensolver context
750: -  i   - iteration number

752:    Output Parameters:
753: .  breakdown - flag indicating that a breakdown has occurred

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

762:    The flag breakdown is set to true if either i=0 and the vector belongs to the
763:    deflation space, or i>0 and the vector is linearly dependent with respect
764:    to the V-vectors.
765: */
766: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,PetscBool *breakdown)
767: {
768:   PetscReal      norm;
769:   PetscBool      lindep;
770:   Vec            w,z;


775:   /* For the first step, use the first initial vector, otherwise a random one */
776:   if (i>0 || eps->nini==0) BVSetRandomColumn(eps->V,i);

778:   /* Force the vector to be in the range of OP for definite generalized problems */
779:   if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
780:     BVCreateVec(eps->V,&w);
781:     BVCopyVec(eps->V,i,w);
782:     BVGetColumn(eps->V,i,&z);
783:     STApply(eps->st,w,z);
784:     BVRestoreColumn(eps->V,i,&z);
785:     VecDestroy(&w);
786:   }

788:   /* Orthonormalize the vector with respect to previous vectors */
789:   BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&lindep);
790:   if (breakdown) *breakdown = lindep;
791:   else if (lindep || norm == 0.0) {
794:   }
795:   BVScaleColumn(eps->V,i,1.0/norm);
796:   PetscFunctionReturn(0);
797: }

799: /*
800:    EPSGetLeftStartVector - Generate a suitable vector to be used as the left starting
801:    vector for the recurrence that builds the left subspace. See EPSGetStartVector().
802: */
803: PetscErrorCode EPSGetLeftStartVector(EPS eps,PetscInt i,PetscBool *breakdown)
804: {
805:   PetscReal      norm;
806:   PetscBool      lindep;


811:   /* For the first step, use the first initial vector, otherwise a random one */
812:   if (i>0 || eps->ninil==0) BVSetRandomColumn(eps->W,i);

814:   /* Orthonormalize the vector with respect to previous vectors */
815:   BVOrthogonalizeColumn(eps->W,i,NULL,&norm,&lindep);
816:   if (breakdown) *breakdown = lindep;
817:   else if (lindep || norm == 0.0) {
819:     SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Unable to generate more left start vectors");
820:   }
821:   BVScaleColumn(eps->W,i,1.0/norm);
822:   PetscFunctionReturn(0);
823: }