Actual source code: nepimpl.h

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

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: #if !defined(_NEPIMPL)
 23: #define _NEPIMPL

 25: #include <slepcnep.h>
 26: #include <slepc-private/slepcimpl.h>

 28: PETSC_EXTERN PetscLogEvent NEP_SetUp,NEP_Solve,NEP_Refine,NEP_FunctionEval,NEP_JacobianEval;

 30: typedef struct _NEPOps *NEPOps;

 32: struct _NEPOps {
 33:   PetscErrorCode (*solve)(NEP);
 34:   PetscErrorCode (*setup)(NEP);
 35:   PetscErrorCode (*setfromoptions)(NEP);
 36:   PetscErrorCode (*publishoptions)(NEP);
 37:   PetscErrorCode (*destroy)(NEP);
 38:   PetscErrorCode (*reset)(NEP);
 39:   PetscErrorCode (*view)(NEP,PetscViewer);
 40:   PetscErrorCode (*computevectors)(NEP);
 41: };

 43: /*
 44:      Maximum number of monitors you can run with a single NEP
 45: */
 46: #define MAXNEPMONITORS 5

 48: typedef enum { NEP_STATE_INITIAL,
 49:                NEP_STATE_SETUP,
 50:                NEP_STATE_SOLVED,
 51:                NEP_STATE_EIGENVECTORS } NEPStateType;

 53: /*
 54:    Defines the NEP data structure.
 55: */
 56: struct _p_NEP {
 57:   PETSCHEADER(struct _NEPOps);
 58:   /*------------------------- User parameters ---------------------------*/
 59:   PetscInt       max_it;           /* maximum number of iterations */
 60:   PetscInt       max_funcs;        /* maximum number of function evaluations */
 61:   PetscInt       nev;              /* number of eigenvalues to compute */
 62:   PetscInt       ncv;              /* number of basis vectors */
 63:   PetscInt       mpd;              /* maximum dimension of projected problem */
 64:   PetscInt       lag;              /* interval to rebuild preconditioner */
 65:   PetscInt       nini;             /* number of initial vectors (negative means not copied yet) */
 66:   PetscScalar    target;           /* target value */
 67:   PetscReal      abstol,rtol,stol; /* user tolerances */
 68:   PetscReal      ktol;             /* tolerance for linear solver */
 69:   PetscBool      cctol;            /* constant correction tolerance */
 70:   PetscReal      ttol;             /* tolerance used in the convergence criterion */
 71:   NEPWhich       which;            /* which part of the spectrum to be sought */
 72:   NEPRefine      refine;           /* type of refinement to be applied after solve */
 73:   PetscReal      reftol;           /* tolerance for refinement */
 74:   PetscInt       rits;             /* number of iterations of the refinement method */
 75:   PetscBool      trackall;         /* whether all the residuals must be computed */

 77:   /*-------------- User-provided functions and contexts -----------------*/
 78:   PetscErrorCode (*computefunction)(NEP,PetscScalar,Mat,Mat,void*);
 79:   PetscErrorCode (*computejacobian)(NEP,PetscScalar,Mat,void*);
 80:   void           *functionctx;
 81:   void           *jacobianctx;
 82:   PetscErrorCode (*converged)(NEP,PetscInt,PetscReal,PetscReal,PetscReal,NEPConvergedReason*,void*);
 83:   PetscErrorCode (*convergeddestroy)(void*);
 84:   void           *convergedctx;
 85:   PetscErrorCode (*monitor[MAXNEPMONITORS])(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt,void*);
 86:   PetscErrorCode (*monitordestroy[MAXNEPMONITORS])(void**);
 87:   void           *monitorcontext[MAXNEPMONITORS];
 88:   PetscInt       numbermonitors;

 90:   /*----------------- Child objects and working data -------------------*/
 91:   DS             ds;               /* direct solver object */
 92:   BV             V;                /* set of basis vectors and computed eigenvectors */
 93:   RG             rg;               /* optional region for filtering */
 94:   PetscRandom    rand;             /* random number generator */
 95:   SlepcSC        sc;               /* sorting criterion data */
 96:   KSP            ksp;              /* linear solver object */
 97:   Mat            function;         /* function matrix */
 98:   Mat            function_pre;     /* function matrix (preconditioner) */
 99:   Mat            jacobian;         /* Jacobian matrix */
100:   Mat            *A;               /* matrix coefficients of split form */
101:   FN             *f;               /* matrix functions of split form */
102:   PetscInt       nt;               /* number of terms in split form */
103:   MatStructure   mstr;             /* pattern of split matrices */
104:   Vec            *IS;              /* references to user-provided initial space */
105:   PetscScalar    *eigr,*eigi;      /* real and imaginary parts of eigenvalues */
106:   PetscReal      *errest;          /* error estimates */
107:   PetscInt       *perm;            /* permutation for eigenvalue ordering */
108:   PetscInt       nwork;            /* number of work vectors */
109:   Vec            *work;            /* work vectors */
110:   void           *data;            /* placeholder for solver-specific stuff */

112:   /* ----------------------- Status variables --------------------------*/
113:   NEPStateType   state;            /* initial -> setup -> solved -> eigenvectors */
114:   PetscInt       nconv;            /* number of converged eigenvalues */
115:   PetscInt       its;              /* number of iterations so far computed */
116:   PetscInt       n,nloc;           /* problem dimensions (global, local) */
117:   PetscInt       nfuncs;           /* number of function evaluations */
118:   PetscBool      split;            /* the nonlinear operator has been set in
119:                                       split form, otherwise user callbacks are used */
120:   NEPConvergedReason reason;
121: };

123: /*
124:     Macros to test valid NEP arguments
125: */
126: #if !defined(PETSC_USE_DEBUG)

128: #define NEPCheckSolved(h,arg) do {} while (0)

130: #else

132: #define NEPCheckSolved(h,arg) \
133:   do { \
134:     if (h->state<NEP_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call NEPSolve() first: Parameter #%d",arg); \
135:   } while (0)

137: #endif

141: PETSC_STATIC_INLINE PetscErrorCode NEP_KSPSolve(NEP nep,Vec b,Vec x)
142: {
144:   PetscInt       lits;

147:   KSPSolve(nep->ksp,b,x);
148:   KSPGetIterationNumber(nep->ksp,&lits);
149:   PetscInfo2(nep,"iter=%D, linear solve iterations=%D\n",nep->its,lits);
150:   return(0);
151: }

153: PETSC_INTERN PetscErrorCode NEPGetDefaultShift(NEP,PetscScalar*);
154: PETSC_INTERN PetscErrorCode NEPComputeResidualNorm_Private(NEP,PetscScalar,Vec,PetscReal*);
155: PETSC_INTERN PetscErrorCode NEPComputeRelativeError_Private(NEP,PetscScalar,Vec,PetscReal*);
156: PETSC_INTERN PetscErrorCode NEPNewtonRefinementSimple(NEP,PetscInt*,PetscReal*,PetscInt);

158: #endif