Actual source code: nepdefault.c

  1: /*
  2:      This file contains some simple default routines for common NEP operations.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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 <slepc-private/nepimpl.h>     /*I "slepcnep.h" I*/
 25: #include <slepcblaslapack.h>

 29: PetscErrorCode NEPReset_Default(NEP nep)
 30: {

 34:   VecDestroyVecs(nep->nwork,&nep->work);
 35:   nep->nwork = 0;
 36:   NEPFreeSolution(nep);
 37:   return(0);
 38: }

 42: /*@
 43:    NEPSetWorkVecs - Sets a number of work vectors into a NEP object

 45:    Collective on NEP

 47:    Input Parameters:
 48: +  nep - nonlinear eigensolver context
 49: -  nw  - number of work vectors to allocate

 51:    Developers Note:
 52:    This is PETSC_EXTERN because it may be required by user plugin NEP
 53:    implementations.

 55:    Level: developer
 56: @*/
 57: PetscErrorCode NEPSetWorkVecs(NEP nep,PetscInt nw)
 58: {

 62:   if (nep->nwork != nw) {
 63:     VecDestroyVecs(nep->nwork,&nep->work);
 64:     nep->nwork = nw;
 65:     VecDuplicateVecs(nep->t,nw,&nep->work);
 66:     PetscLogObjectParents(nep,nw,nep->work);
 67:   }
 68:   return(0);
 69: }

 73: /*
 74:   NEPGetDefaultShift - Return the value of sigma to start the nonlinear iteration.
 75:  */
 76: PetscErrorCode NEPGetDefaultShift(NEP nep,PetscScalar *sigma)
 77: {
 80:   switch (nep->which) {
 81:     case NEP_LARGEST_MAGNITUDE:
 82:     case NEP_LARGEST_IMAGINARY:
 83:       *sigma = 1.0;   /* arbitrary value */
 84:       break;
 85:     case NEP_SMALLEST_MAGNITUDE:
 86:     case NEP_SMALLEST_IMAGINARY:
 87:       *sigma = 0.0;
 88:       break;
 89:     case NEP_LARGEST_REAL:
 90:       *sigma = PETSC_MAX_REAL;
 91:       break;
 92:     case NEP_SMALLEST_REAL:
 93:       *sigma = PETSC_MIN_REAL;
 94:       break;
 95:     case NEP_TARGET_MAGNITUDE:
 96:     case NEP_TARGET_REAL:
 97:     case NEP_TARGET_IMAGINARY:
 98:       *sigma = nep->target;
 99:       break;
100:   }
101:   return(0);
102: }

106: /*
107:   NEPConvergedDefault - Checks convergence of the nonlinear eigensolver.
108: */
109: PetscErrorCode NEPConvergedDefault(NEP nep,PetscInt it,PetscReal xnorm,PetscReal snorm,PetscReal fnorm,NEPConvergedReason *reason,void *ctx)
110: {


117:   *reason = NEP_CONVERGED_ITERATING;

119:   if (!it) {
120:     /* set parameter for default relative tolerance convergence test */
121:     nep->ttol = fnorm*nep->rtol;
122:   }
123:   if (PetscIsInfOrNanReal(fnorm)) {
124:     PetscInfo(nep,"Failed to converged, function norm is NaN\n");
125:     *reason = NEP_DIVERGED_FNORM_NAN;
126:   } else if (fnorm < nep->abstol) {
127:     PetscInfo2(nep,"Converged due to function norm %14.12e < %14.12e\n",(double)fnorm,(double)nep->abstol);
128:     *reason = NEP_CONVERGED_FNORM_ABS;
129:   } else if (nep->nfuncs >= nep->max_funcs) {
130:     PetscInfo2(nep,"Exceeded maximum number of function evaluations: %D > %D\n",nep->nfuncs,nep->max_funcs);
131:     *reason = NEP_DIVERGED_FUNCTION_COUNT;
132:   }

134:   if (it && !*reason) {
135:     if (fnorm <= nep->ttol) {
136:       PetscInfo2(nep,"Converged due to function norm %14.12e < %14.12e (relative tolerance)\n",(double)fnorm,(double)nep->ttol);
137:       *reason = NEP_CONVERGED_FNORM_RELATIVE;
138:     } else if (snorm < nep->stol*xnorm) {
139:       PetscInfo3(nep,"Converged due to small update length: %14.12e < %14.12e * %14.12e\n",(double)snorm,(double)nep->stol,(double)xnorm);
140:       *reason = NEP_CONVERGED_SNORM_RELATIVE;
141:     }
142:   }
143:   return(0);
144: }