Actual source code: nepdefault.c

slepc-3.18.0 2022-10-01
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:    Simple default routines for common NEP operations
 12: */

 14: #include <slepc/private/nepimpl.h>

 16: /*@
 17:    NEPSetWorkVecs - Sets a number of work vectors into a NEP object

 19:    Collective on nep

 21:    Input Parameters:
 22: +  nep - nonlinear eigensolver context
 23: -  nw  - number of work vectors to allocate

 25:    Developer Notes:
 26:    This is SLEPC_EXTERN because it may be required by user plugin NEP
 27:    implementations.

 29:    Level: developer

 31: .seealso: NEPSetUp()
 32: @*/
 33: PetscErrorCode NEPSetWorkVecs(NEP nep,PetscInt nw)
 34: {
 35:   Vec            t;

 37:   if (nep->nwork < nw) {
 38:     VecDestroyVecs(nep->nwork,&nep->work);
 39:     nep->nwork = nw;
 40:     BVGetColumn(nep->V,0,&t);
 41:     VecDuplicateVecs(t,nw,&nep->work);
 42:     BVRestoreColumn(nep->V,0,&t);
 43:   }
 44:   return 0;
 45: }

 47: /*
 48:   NEPGetDefaultShift - Return the value of sigma to start the nonlinear iteration.
 49:  */
 50: PetscErrorCode NEPGetDefaultShift(NEP nep,PetscScalar *sigma)
 51: {
 53:   switch (nep->which) {
 54:     case NEP_LARGEST_MAGNITUDE:
 55:     case NEP_LARGEST_IMAGINARY:
 56:     case NEP_ALL:
 57:     case NEP_WHICH_USER:
 58:       *sigma = 1.0;   /* arbitrary value */
 59:       break;
 60:     case NEP_SMALLEST_MAGNITUDE:
 61:     case NEP_SMALLEST_IMAGINARY:
 62:       *sigma = 0.0;
 63:       break;
 64:     case NEP_LARGEST_REAL:
 65:       *sigma = PETSC_MAX_REAL;
 66:       break;
 67:     case NEP_SMALLEST_REAL:
 68:       *sigma = PETSC_MIN_REAL;
 69:       break;
 70:     case NEP_TARGET_MAGNITUDE:
 71:     case NEP_TARGET_REAL:
 72:     case NEP_TARGET_IMAGINARY:
 73:       *sigma = nep->target;
 74:       break;
 75:   }
 76:   return 0;
 77: }

 79: /*
 80:   NEPConvergedRelative - Checks convergence relative to the eigenvalue.
 81: */
 82: PetscErrorCode NEPConvergedRelative(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 83: {
 84:   PetscReal w;

 86:   w = SlepcAbsEigenvalue(eigr,eigi);
 87:   *errest = res/w;
 88:   return 0;
 89: }

 91: /*
 92:   NEPConvergedAbsolute - Checks convergence absolutely.
 93: */
 94: PetscErrorCode NEPConvergedAbsolute(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 95: {
 96:   *errest = res;
 97:   return 0;
 98: }

100: /*
101:   NEPConvergedNorm - Checks convergence relative to the matrix norms.
102: */
103: PetscErrorCode NEPConvergedNorm(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
104: {
105:   PetscScalar    s;
106:   PetscReal      w=0.0;
107:   PetscInt       j;
108:   PetscBool      flg;

110:   if (nep->fui!=NEP_USER_INTERFACE_SPLIT) {
111:     NEPComputeFunction(nep,eigr,nep->function,nep->function);
112:     MatHasOperation(nep->function,MATOP_NORM,&flg);
114:     MatNorm(nep->function,NORM_INFINITY,&w);
115:   } else {
116:     /* initialization of matrix norms */
117:     if (!nep->nrma[0]) {
118:       for (j=0;j<nep->nt;j++) {
119:         MatHasOperation(nep->A[j],MATOP_NORM,&flg);
121:         MatNorm(nep->A[j],NORM_INFINITY,&nep->nrma[j]);
122:       }
123:     }
124:     for (j=0;j<nep->nt;j++) {
125:       FNEvaluateFunction(nep->f[j],eigr,&s);
126:       w = w + nep->nrma[j]*PetscAbsScalar(s);
127:     }
128:   }
129:   *errest = res/w;
130:   return 0;
131: }

133: /*@C
134:    NEPStoppingBasic - Default routine to determine whether the outer eigensolver
135:    iteration must be stopped.

137:    Collective on nep

139:    Input Parameters:
140: +  nep    - nonlinear eigensolver context obtained from NEPCreate()
141: .  its    - current number of iterations
142: .  max_it - maximum number of iterations
143: .  nconv  - number of currently converged eigenpairs
144: .  nev    - number of requested eigenpairs
145: -  ctx    - context (not used here)

147:    Output Parameter:
148: .  reason - result of the stopping test

150:    Notes:
151:    A positive value of reason indicates that the iteration has finished successfully
152:    (converged), and a negative value indicates an error condition (diverged). If
153:    the iteration needs to be continued, reason must be set to NEP_CONVERGED_ITERATING
154:    (zero).

156:    NEPStoppingBasic() will stop if all requested eigenvalues are converged, or if
157:    the maximum number of iterations has been reached.

159:    Use NEPSetStoppingTest() to provide your own test instead of using this one.

161:    Level: advanced

163: .seealso: NEPSetStoppingTest(), NEPConvergedReason, NEPGetConvergedReason()
164: @*/
165: PetscErrorCode NEPStoppingBasic(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx)
166: {
167:   *reason = NEP_CONVERGED_ITERATING;
168:   if (nconv >= nev) {
169:     PetscInfo(nep,"Nonlinear eigensolver finished successfully: %" PetscInt_FMT " eigenpairs converged at iteration %" PetscInt_FMT "\n",nconv,its);
170:     *reason = NEP_CONVERGED_TOL;
171:   } else if (its >= max_it) {
172:     *reason = NEP_DIVERGED_ITS;
173:     PetscInfo(nep,"Nonlinear eigensolver iteration reached maximum number of iterations (%" PetscInt_FMT ")\n",its);
174:   }
175:   return 0;
176: }

178: PetscErrorCode NEPComputeVectors_Schur(NEP nep)
179: {
180:   Mat            Z;

182:   DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
183:   DSGetMat(nep->ds,DS_MAT_X,&Z);
184:   BVMultInPlace(nep->V,Z,0,nep->nconv);
185:   DSRestoreMat(nep->ds,DS_MAT_X,&Z);
186:   BVNormalize(nep->V,nep->eigi);
187:   return 0;
188: }