Actual source code: nepdefault.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:    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:     PetscLogObjectParents(nep,nw,nep->work);
 44:   }
 45:   PetscFunctionReturn(0);
 46: }

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

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

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

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

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

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

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

138:    Collective on nep

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

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

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

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

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

162:    Level: advanced

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

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

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