Actual source code: epssolve.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: {

 23:   EPSCheckSolved(eps,1);
 24:   if (eps->state==EPS_STATE_SOLVED && eps->ops->computevectors) {
 25:     (*eps->ops->computevectors)(eps);
 26:   }
 27:   eps->state = EPS_STATE_EIGENVECTORS;
 28:   return(0);
 29: }

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

 33: static PetscErrorCode EPSComputeValues(EPS eps)
 34: {
 36:   PetscBool      injective,iscomp,isfilter;
 37:   PetscInt       i,n,aux,nconv0;
 38:   Mat            A,B=NULL,G,Z;

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

100: /*@
101:    EPSSolve - Solves the eigensystem.

103:    Collective on eps

105:    Input Parameter:
106: .  eps - eigensolver context obtained from EPSCreate()

108:    Options Database Keys:
109: +  -eps_view - print information about the solver used
110: .  -eps_view_mat0 binary - save the first matrix (A) to the default binary viewer
111: .  -eps_view_mat1 binary - save the second matrix (B) to the default binary viewer
112: .  -eps_view_vectors binary - save the computed eigenvectors to the default binary viewer
113: .  -eps_view_values - print computed eigenvalues
114: .  -eps_converged_reason - print reason for convergence, and number of iterations
115: .  -eps_error_absolute - print absolute errors of each eigenpair
116: .  -eps_error_relative - print relative errors of each eigenpair
117: -  -eps_error_backward - print backward errors of each eigenpair

119:    Level: beginner

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

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

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

148:   /* Call solver */
149:   (*eps->ops->solve)(eps);
150:   if (!eps->reason) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
151:   eps->state = EPS_STATE_SOLVED;

153:   /* Only the first nconv columns contain useful information (except in CISS) */
154:   BVSetActiveColumns(eps->V,0,eps->nconv);
155:   if (eps->twosided) { BVSetActiveColumns(eps->W,0,eps->nconv); }

157:   /* If inplace, purify eigenvectors before reverting operator */
158:   STGetMatMode(eps->st,&matmode);
159:   if (matmode == ST_MATMODE_INPLACE && eps->ispositive) {
160:     EPSComputeVectors(eps);
161:   }
162:   STPostSolve(eps->st);

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

167: #if !defined(PETSC_USE_COMPLEX)
168:   /* Reorder conjugate eigenvalues (positive imaginary first) */
169:   for (i=0;i<eps->nconv-1;i++) {
170:     if (eps->eigi[i] != 0) {
171:       if (eps->eigi[i] < 0) {
172:         eps->eigi[i] = -eps->eigi[i];
173:         eps->eigi[i+1] = -eps->eigi[i+1];
174:         /* the next correction only works with eigenvectors */
175:         EPSComputeVectors(eps);
176:         BVScaleColumn(eps->V,i+1,-1.0);
177:       }
178:       i++;
179:     }
180:   }
181: #endif

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

187:   /* Various viewers */
188:   EPSViewFromOptions(eps,NULL,"-eps_view");
189:   EPSConvergedReasonViewFromOptions(eps);
190:   EPSErrorViewFromOptions(eps);
191:   EPSValuesViewFromOptions(eps);
192:   EPSVectorsViewFromOptions(eps);
193:   EPSGetOperators(eps,&A,&B);
194:   MatViewFromOptions(A,(PetscObject)eps,"-eps_view_mat0");
195:   if (eps->isgeneralized) {
196:     MatViewFromOptions(B,(PetscObject)eps,"-eps_view_mat1");
197:   }

199:   /* Remove deflation and initial subspaces */
200:   if (eps->nds) {
201:     BVSetNumConstraints(eps->V,0);
202:     eps->nds = 0;
203:   }
204:   eps->nini = 0;
205:   return(0);
206: }

208: /*@
209:    EPSGetIterationNumber - Gets the current iteration number. If the
210:    call to EPSSolve() is complete, then it returns the number of iterations
211:    carried out by the solution method.

213:    Not Collective

215:    Input Parameter:
216: .  eps - the eigensolver context

218:    Output Parameter:
219: .  its - number of iterations

221:    Note:
222:    During the i-th iteration this call returns i-1. If EPSSolve() is
223:    complete, then parameter "its" contains either the iteration number at
224:    which convergence was successfully reached, or failure was detected.
225:    Call EPSGetConvergedReason() to determine if the solver converged or
226:    failed and why.

228:    Level: intermediate

230: .seealso: EPSGetConvergedReason(), EPSSetTolerances()
231: @*/
232: PetscErrorCode EPSGetIterationNumber(EPS eps,PetscInt *its)
233: {
237:   *its = eps->its;
238:   return(0);
239: }

241: /*@
242:    EPSGetConverged - Gets the number of converged eigenpairs.

244:    Not Collective

246:    Input Parameter:
247: .  eps - the eigensolver context

249:    Output Parameter:
250: .  nconv - number of converged eigenpairs

252:    Note:
253:    This function should be called after EPSSolve() has finished.

255:    Level: beginner

257: .seealso: EPSSetDimensions(), EPSSolve()
258: @*/
259: PetscErrorCode EPSGetConverged(EPS eps,PetscInt *nconv)
260: {
264:   EPSCheckSolved(eps,1);
265:   *nconv = eps->nconv;
266:   return(0);
267: }

269: /*@
270:    EPSGetConvergedReason - Gets the reason why the EPSSolve() iteration was
271:    stopped.

273:    Not Collective

275:    Input Parameter:
276: .  eps - the eigensolver context

278:    Output Parameter:
279: .  reason - negative value indicates diverged, positive value converged

281:    Notes:
282:    Possible values for reason are
283: +  EPS_CONVERGED_TOL - converged up to tolerance
284: .  EPS_CONVERGED_USER - converged due to a user-defined condition
285: .  EPS_DIVERGED_ITS - required more than max_it iterations to reach convergence
286: .  EPS_DIVERGED_BREAKDOWN - generic breakdown in method
287: -  EPS_DIVERGED_SYMMETRY_LOST - pseudo-Lanczos was not able to keep symmetry

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

291:    Level: intermediate

293: .seealso: EPSSetTolerances(), EPSSolve(), EPSConvergedReason
294: @*/
295: PetscErrorCode EPSGetConvergedReason(EPS eps,EPSConvergedReason *reason)
296: {
300:   EPSCheckSolved(eps,1);
301:   *reason = eps->reason;
302:   return(0);
303: }

305: /*@
306:    EPSGetInvariantSubspace - Gets an orthonormal basis of the computed invariant
307:    subspace.

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

311:    Input Parameter:
312: .  eps - the eigensolver context

314:    Output Parameter:
315: .  v - an array of vectors

317:    Notes:
318:    This function should be called after EPSSolve() has finished.

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

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

330:    Level: intermediate

332: .seealso: EPSGetEigenpair(), EPSGetConverged(), EPSSolve()
333: @*/
334: PetscErrorCode EPSGetInvariantSubspace(EPS eps,Vec v[])
335: {
337:   PetscInt       i;
338:   BV             V=eps->V;
339:   Vec            w;

345:   EPSCheckSolved(eps,1);
346:   if (!eps->ishermitian && eps->state==EPS_STATE_EIGENVECTORS) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSGetInvariantSubspace must be called before EPSGetEigenpair,EPSGetEigenvector or EPSComputeError");
347:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
348:     BVDuplicateResize(eps->V,eps->nconv,&V);
349:     BVSetActiveColumns(eps->V,0,eps->nconv);
350:     BVCopy(eps->V,V);
351:     for (i=0;i<eps->nconv;i++) {
352:       BVGetColumn(V,i,&w);
353:       VecPointwiseDivide(w,w,eps->D);
354:       BVRestoreColumn(V,i,&w);
355:     }
356:     BVOrthogonalize(V,NULL);
357:   }
358:   for (i=0;i<eps->nconv;i++) {
359:     BVCopyVec(V,i,v[i]);
360:   }
361:   if (eps->balance!=EPS_BALANCE_NONE && eps->D) {
362:     BVDestroy(&V);
363:   }
364:   return(0);
365: }

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

371:    Logically Collective on eps

373:    Input Parameters:
374: +  eps - eigensolver context
375: -  i   - index of the solution

377:    Output Parameters:
378: +  eigr - real part of eigenvalue
379: .  eigi - imaginary part of eigenvalue
380: .  Vr   - real part of eigenvector
381: -  Vi   - imaginary part of eigenvector

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

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

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

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

401:    Level: beginner

403: .seealso: EPSGetEigenvalue(), EPSGetEigenvector(), EPSGetLeftEigenvector(), EPSSolve(),
404:           EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetInvariantSubspace()
405: @*/
406: PetscErrorCode EPSGetEigenpair(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)
407: {

413:   EPSCheckSolved(eps,1);
414:   if (i<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
415:   if (i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()");
416:   EPSGetEigenvalue(eps,i,eigr,eigi);
417:   if (Vr || Vi) { EPSGetEigenvector(eps,i,Vr,Vi); }
418:   return(0);
419: }

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

424:    Not Collective

426:    Input Parameters:
427: +  eps - eigensolver context
428: -  i   - index of the solution

430:    Output Parameters:
431: +  eigr - real part of eigenvalue
432: -  eigi - imaginary part of eigenvalue

434:    Notes:
435:    If the eigenvalue is real, then eigi is set to zero. If PETSc is
436:    configured with complex scalars the eigenvalue is stored
437:    directly in eigr (eigi is set to zero).

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

443:    Level: beginner

445: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair()
446: @*/
447: PetscErrorCode EPSGetEigenvalue(EPS eps,PetscInt i,PetscScalar *eigr,PetscScalar *eigi)
448: {
449:   PetscInt k;

453:   EPSCheckSolved(eps,1);
454:   if (i<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
455:   if (i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()");
456:   k = eps->perm[i];
457: #if defined(PETSC_USE_COMPLEX)
458:   if (eigr) *eigr = eps->eigr[k];
459:   if (eigi) *eigi = 0;
460: #else
461:   if (eigr) *eigr = eps->eigr[k];
462:   if (eigi) *eigi = eps->eigi[k];
463: #endif
464:   return(0);
465: }

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

470:    Logically Collective on eps

472:    Input Parameters:
473: +  eps - eigensolver context
474: -  i   - index of the solution

476:    Output Parameters:
477: +  Vr   - real part of eigenvector
478: -  Vi   - imaginary part of eigenvector

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

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

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

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

497:    Level: beginner

499: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSGetEigenpair(), EPSGetLeftEigenvector()
500: @*/
501: PetscErrorCode EPSGetEigenvector(EPS eps,PetscInt i,Vec Vr,Vec Vi)
502: {
504:   PetscInt       k;

511:   EPSCheckSolved(eps,1);
512:   if (i<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index cannot be negative");
513:   if (i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The index can be nconv-1 at most, see EPSGetConverged()");
514:   EPSComputeVectors(eps);
515:   k = eps->perm[i];
516:   BV_GetEigenvector(eps->V,k,eps->eigi[k],Vr,Vi);
517:   return(0);
518: }

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

523:    Logically Collective on eps

525:    Input Parameters:
526: +  eps - eigensolver context
527: -  i   - index of the solution

529:    Output Parameters:
530: +  Wr   - real part of left eigenvector
531: -  Wi   - imaginary part of left eigenvector

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

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

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

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

549:    Level: intermediate

551: .seealso: EPSGetEigenvector(), EPSGetConverged(), EPSSetWhichEigenpairs(), EPSSetTwoSided()
552: @*/
553: PetscErrorCode EPSGetLeftEigenvector(EPS eps,PetscInt i,Vec Wr,Vec Wi)
554: {
556:   PetscInt       k;

563:   EPSCheckSolved(eps,1);
564:   if (!eps->twosided) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"Must request left vectors with EPSSetTwoSided");
565:   if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
566:   EPSComputeVectors(eps);
567:   k = eps->perm[i];
568:   BV_GetEigenvector(eps->W,k,eps->eigi[k],Wr,Wi);
569:   return(0);
570: }

572: /*@
573:    EPSGetErrorEstimate - Returns the error estimate associated to the i-th
574:    computed eigenpair.

576:    Not Collective

578:    Input Parameter:
579: +  eps - eigensolver context
580: -  i   - index of eigenpair

582:    Output Parameter:
583: .  errest - the error estimate

585:    Notes:
586:    This is the error estimate used internally by the eigensolver. The actual
587:    error bound can be computed with EPSComputeError(). See also the users
588:    manual for details.

590:    Level: advanced

592: .seealso: EPSComputeError()
593: @*/
594: PetscErrorCode EPSGetErrorEstimate(EPS eps,PetscInt i,PetscReal *errest)
595: {
599:   EPSCheckSolved(eps,1);
600:   if (i<0 || i>=eps->nconv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
601:   *errest = eps->errest[eps->perm[i]];
602:   return(0);
603: }

605: /*
606:    EPSComputeResidualNorm_Private - Computes the norm of the residual vector
607:    associated with an eigenpair.

609:    Input Parameters:
610:      trans - whether A' must be used instead of A
611:      kr,ki - eigenvalue
612:      xr,xi - eigenvector
613:      z     - three work vectors (the second one not referenced in complex scalars)
614: */
615: PetscErrorCode EPSComputeResidualNorm_Private(EPS eps,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec *z,PetscReal *norm)
616: {
618:   PetscInt       nmat;
619:   Mat            A,B;
620:   Vec            u,w;
621:   PetscScalar    alpha;
622: #if !defined(PETSC_USE_COMPLEX)
623:   Vec            v;
624:   PetscReal      ni,nr;
625: #endif
626:   PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultHermitianTranspose: MatMult;

629:   u = z[0]; w = z[2];
630:   STGetNumMatrices(eps->st,&nmat);
631:   STGetMatrix(eps->st,0,&A);
632:   if (nmat>1) { STGetMatrix(eps->st,1,&B); }

634: #if !defined(PETSC_USE_COMPLEX)
635:   v = z[1];
636:   if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
637: #endif
638:     (*matmult)(A,xr,u);                          /* u=A*x */
639:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
640:       if (nmat>1) { (*matmult)(B,xr,w); }
641:       else { VecCopy(xr,w); }                    /* w=B*x */
642:       alpha = trans? -PetscConj(kr): -kr;
643:       VecAXPY(u,alpha,w);                        /* u=A*x-k*B*x */
644:     }
645:     VecNorm(u,NORM_2,norm);
646: #if !defined(PETSC_USE_COMPLEX)
647:   } else {
648:     (*matmult)(A,xr,u);                          /* u=A*xr */
649:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
650:       if (nmat>1) { (*matmult)(B,xr,v); }
651:       else { VecCopy(xr,v); }                    /* v=B*xr */
652:       VecAXPY(u,-kr,v);                          /* u=A*xr-kr*B*xr */
653:       if (nmat>1) { (*matmult)(B,xi,w); }
654:       else { VecCopy(xi,w); }                    /* w=B*xi */
655:       VecAXPY(u,trans?-ki:ki,w);                 /* u=A*xr-kr*B*xr+ki*B*xi */
656:     }
657:     VecNorm(u,NORM_2,&nr);
658:     (*matmult)(A,xi,u);                          /* u=A*xi */
659:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
660:       VecAXPY(u,-kr,w);                          /* u=A*xi-kr*B*xi */
661:       VecAXPY(u,trans?ki:-ki,v);                 /* u=A*xi-kr*B*xi-ki*B*xr */
662:     }
663:     VecNorm(u,NORM_2,&ni);
664:     *norm = SlepcAbsEigenvalue(nr,ni);
665:   }
666: #endif
667:   return(0);
668: }

670: /*@
671:    EPSComputeError - Computes the error (based on the residual norm) associated
672:    with the i-th computed eigenpair.

674:    Collective on eps

676:    Input Parameter:
677: +  eps  - the eigensolver context
678: .  i    - the solution index
679: -  type - the type of error to compute

681:    Output Parameter:
682: .  error - the error

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

688:    Level: beginner

690: .seealso: EPSErrorType, EPSSolve(), EPSGetErrorEstimate()
691: @*/
692: PetscErrorCode EPSComputeError(EPS eps,PetscInt i,EPSErrorType type,PetscReal *error)
693: {
695:   Mat            A,B;
696:   Vec            xr,xi,w[3];
697:   PetscReal      t,vecnorm=1.0,errorl;
698:   PetscScalar    kr,ki;
699:   PetscBool      flg;

706:   EPSCheckSolved(eps,1);

708:   /* allocate work vectors */
709: #if defined(PETSC_USE_COMPLEX)
710:   EPSSetWorkVecs(eps,3);
711:   xi   = NULL;
712:   w[1] = NULL;
713: #else
714:   EPSSetWorkVecs(eps,5);
715:   xi   = eps->work[3];
716:   w[1] = eps->work[4];
717: #endif
718:   xr   = eps->work[0];
719:   w[0] = eps->work[1];
720:   w[2] = eps->work[2];

722:   /* compute residual norm */
723:   EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
724:   EPSComputeResidualNorm_Private(eps,PETSC_FALSE,kr,ki,xr,xi,w,error);

726:   /* compute 2-norm of eigenvector */
727:   if (eps->problem_type==EPS_GHEP) {
728:     VecNorm(xr,NORM_2,&vecnorm);
729:   }

731:   /* if two-sided, compute left residual norm and take the maximum */
732:   if (eps->twosided) {
733:     EPSGetLeftEigenvector(eps,i,xr,xi);
734:     EPSComputeResidualNorm_Private(eps,PETSC_TRUE,kr,ki,xr,xi,w,&errorl);
735:     *error = PetscMax(*error,errorl);
736:   }

738:   /* compute error */
739:   switch (type) {
740:     case EPS_ERROR_ABSOLUTE:
741:       break;
742:     case EPS_ERROR_RELATIVE:
743:       *error /= SlepcAbsEigenvalue(kr,ki)*vecnorm;
744:       break;
745:     case EPS_ERROR_BACKWARD:
746:       /* initialization of matrix norms */
747:       if (!eps->nrma) {
748:         STGetMatrix(eps->st,0,&A);
749:         MatHasOperation(A,MATOP_NORM,&flg);
750:         if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
751:         MatNorm(A,NORM_INFINITY,&eps->nrma);
752:       }
753:       if (eps->isgeneralized) {
754:         if (!eps->nrmb) {
755:           STGetMatrix(eps->st,1,&B);
756:           MatHasOperation(B,MATOP_NORM,&flg);
757:           if (!flg) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
758:           MatNorm(B,NORM_INFINITY,&eps->nrmb);
759:         }
760:       } else eps->nrmb = 1.0;
761:       t = SlepcAbsEigenvalue(kr,ki);
762:       *error /= (eps->nrma+t*eps->nrmb)*vecnorm;
763:       break;
764:     default:
765:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
766:   }
767:   return(0);
768: }

770: /*
771:    EPSGetStartVector - Generate a suitable vector to be used as the starting vector
772:    for the recurrence that builds the right subspace.

774:    Collective on eps

776:    Input Parameters:
777: +  eps - the eigensolver context
778: -  i   - iteration number

780:    Output Parameters:
781: .  breakdown - flag indicating that a breakdown has occurred

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

790:    The flag breakdown is set to true if either i=0 and the vector belongs to the
791:    deflation space, or i>0 and the vector is linearly dependent with respect
792:    to the V-vectors.
793: */
794: PetscErrorCode EPSGetStartVector(EPS eps,PetscInt i,PetscBool *breakdown)
795: {
797:   PetscReal      norm;
798:   PetscBool      lindep;
799:   Vec            w,z;


805:   /* For the first step, use the first initial vector, otherwise a random one */
806:   if (i>0 || eps->nini==0) {
807:     BVSetRandomColumn(eps->V,i);
808:   }

810:   /* Force the vector to be in the range of OP for definite generalized problems */
811:   if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
812:     BVCreateVec(eps->V,&w);
813:     BVCopyVec(eps->V,i,w);
814:     BVGetColumn(eps->V,i,&z);
815:     STApply(eps->st,w,z);
816:     BVRestoreColumn(eps->V,i,&z);
817:     VecDestroy(&w);
818:   }

820:   /* Orthonormalize the vector with respect to previous vectors */
821:   BVOrthogonalizeColumn(eps->V,i,NULL,&norm,&lindep);
822:   if (breakdown) *breakdown = lindep;
823:   else if (lindep || norm == 0.0) {
824:     if (i==0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Initial vector is zero or belongs to the deflation space");
825:     else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to generate more start vectors");
826:   }
827:   BVScaleColumn(eps->V,i,1.0/norm);
828:   return(0);
829: }

831: /*
832:    EPSGetLeftStartVector - Generate a suitable vector to be used as the left starting
833:    vector for the recurrence that builds the left subspace. See EPSGetStartVector().
834: */
835: PetscErrorCode EPSGetLeftStartVector(EPS eps,PetscInt i,PetscBool *breakdown)
836: {
838:   PetscReal      norm;
839:   PetscBool      lindep;


845:   /* For the first step, use the first initial vector, otherwise a random one */
846:   if (i>0 || eps->ninil==0) {
847:     BVSetRandomColumn(eps->W,i);
848:   }

850:   /* Orthonormalize the vector with respect to previous vectors */
851:   BVOrthogonalizeColumn(eps->W,i,NULL,&norm,&lindep);
852:   if (breakdown) *breakdown = lindep;
853:   else if (lindep || norm == 0.0) {
854:     if (i==0) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Left initial vector is zero");
855:     else SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unable to generate more left start vectors");
856:   }
857:   BVScaleColumn(eps->W,i,1.0/norm);
858:   return(0);
859: }