Actual source code: epsimpl.h

slepc-3.7.4 2017-05-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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(_EPSIMPL)
 23: #define _EPSIMPL

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

 28: PETSC_EXTERN PetscBool EPSRegisterAllCalled;
 29: PETSC_EXTERN PetscErrorCode EPSRegisterAll(void);
 30: PETSC_EXTERN PetscLogEvent EPS_SetUp,EPS_Solve;

 32: typedef struct _EPSOps *EPSOps;

 34: struct _EPSOps {
 35:   PetscErrorCode (*solve)(EPS);
 36:   PetscErrorCode (*setup)(EPS);
 37:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,EPS);
 38:   PetscErrorCode (*publishoptions)(EPS);
 39:   PetscErrorCode (*destroy)(EPS);
 40:   PetscErrorCode (*reset)(EPS);
 41:   PetscErrorCode (*view)(EPS,PetscViewer);
 42:   PetscErrorCode (*backtransform)(EPS);
 43:   PetscErrorCode (*computevectors)(EPS);
 44: };

 46: /*
 47:      Maximum number of monitors you can run with a single EPS
 48: */
 49: #define MAXEPSMONITORS 5

 51: typedef enum { EPS_STATE_INITIAL,
 52:                EPS_STATE_SETUP,
 53:                EPS_STATE_SOLVED,
 54:                EPS_STATE_EIGENVECTORS } EPSStateType;

 56: /*
 57:    Defines the EPS data structure.
 58: */
 59: struct _p_EPS {
 60:   PETSCHEADER(struct _EPSOps);
 61:   /*------------------------- User parameters ---------------------------*/
 62:   PetscInt       max_it;           /* maximum number of iterations */
 63:   PetscInt       nev;              /* number of eigenvalues to compute */
 64:   PetscInt       ncv;              /* number of basis vectors */
 65:   PetscInt       mpd;              /* maximum dimension of projected problem */
 66:   PetscInt       nini;             /* number of initial vectors (negative means not copied yet) */
 67:   PetscInt       nds;              /* number of basis vectors of deflation space */
 68:   PetscScalar    target;           /* target value */
 69:   PetscReal      tol;              /* tolerance */
 70:   EPSConv        conv;             /* convergence test */
 71:   EPSStop        stop;             /* stopping test */
 72:   EPSWhich       which;            /* which part of the spectrum to be sought */
 73:   PetscReal      inta,intb;        /* interval [a,b] for spectrum slicing */
 74:   EPSProblemType problem_type;     /* which kind of problem to be solved */
 75:   EPSExtraction  extraction;       /* which kind of extraction to be applied */
 76:   EPSBalance     balance;          /* the balancing method */
 77:   PetscInt       balance_its;      /* number of iterations of the balancing method */
 78:   PetscReal      balance_cutoff;   /* cutoff value for balancing */
 79:   PetscBool      trueres;          /* whether the true residual norm must be computed */
 80:   PetscBool      trackall;         /* whether all the residuals must be computed */
 81:   PetscBool      purify;           /* whether eigenvectors need to be purified */

 83:   /*-------------- User-provided functions and contexts -----------------*/
 84:   PetscErrorCode (*converged)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 85:   PetscErrorCode (*convergeddestroy)(void*);
 86:   PetscErrorCode (*stopping)(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
 87:   PetscErrorCode (*stoppingdestroy)(void*);
 88:   PetscErrorCode (*arbitrary)(PetscScalar,PetscScalar,Vec,Vec,PetscScalar*,PetscScalar*,void*);
 89:   void           *convergedctx;
 90:   void           *stoppingctx;
 91:   void           *arbitraryctx;
 92:   PetscErrorCode (*monitor[MAXEPSMONITORS])(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
 93:   PetscErrorCode (*monitordestroy[MAXEPSMONITORS])(void**);
 94:   void           *monitorcontext[MAXEPSMONITORS];
 95:   PetscInt       numbermonitors;

 97:   /*----------------- Child objects and working data -------------------*/
 98:   ST             st;               /* spectral transformation object */
 99:   DS             ds;               /* direct solver object */
100:   BV             V;                /* set of basis vectors and computed eigenvectors */
101:   RG             rg;               /* optional region for filtering */
102:   SlepcSC        sc;               /* sorting criterion data */
103:   Vec            D;                /* diagonal matrix for balancing */
104:   Vec            *IS;              /* references to user-provided initial space */
105:   Vec            *defl;            /* references to user-provided deflation space */
106:   PetscScalar    *eigr,*eigi;      /* real and imaginary parts of eigenvalues */
107:   PetscReal      *errest;          /* error estimates */
108:   PetscScalar    *rr,*ri;          /* values computed by user's arbitrary selection function */
109:   PetscInt       *perm;            /* permutation for eigenvalue ordering */
110:   PetscInt       nwork;            /* number of work vectors */
111:   Vec            *work;            /* work vectors */
112:   void           *data;            /* placeholder for solver-specific stuff */

114:   /* ----------------------- Status variables --------------------------*/
115:   EPSStateType   state;            /* initial -> setup -> solved -> eigenvectors */
116:   PetscInt       nconv;            /* number of converged eigenvalues */
117:   PetscInt       its;              /* number of iterations so far computed */
118:   PetscInt       n,nloc;           /* problem dimensions (global, local) */
119:   PetscReal      nrma,nrmb;        /* computed matrix norms */
120:   PetscBool      isgeneralized;
121:   PetscBool      ispositive;
122:   PetscBool      ishermitian;
123:   EPSConvergedReason reason;
124: };

126: /*
127:     Macros to test valid EPS arguments
128: */
129: #if !defined(PETSC_USE_DEBUG)

131: #define EPSCheckSolved(h,arg) do {} while (0)

133: #else

135: #define EPSCheckSolved(h,arg) \
136:   do { \
137:     if (h->state<EPS_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call EPSSolve() first: Parameter #%d",arg); \
138:   } while (0)

140: #endif

144: /*
145:   EPS_SetInnerProduct - set B matrix for inner product if appropriate.
146: */
147: PETSC_STATIC_INLINE PetscErrorCode EPS_SetInnerProduct(EPS eps)
148: {
150:   Mat            B;

153:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
154:   if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
155:     STGetBilinearForm(eps->st,&B);
156:     BVSetMatrix(eps->V,B,PetscNot(eps->ispositive));
157:     MatDestroy(&B);
158:   } else {
159:     BVSetMatrix(eps->V,NULL,PETSC_FALSE);
160:   }
161:   return(0);
162: }

164: PETSC_INTERN PetscErrorCode EPSSetWhichEigenpairs_Default(EPS);
165: PETSC_INTERN PetscErrorCode EPSSetDimensions_Default(EPS,PetscInt,PetscInt*,PetscInt*);
166: PETSC_INTERN PetscErrorCode EPSBackTransform_Default(EPS);
167: PETSC_INTERN PetscErrorCode EPSComputeVectors(EPS);
168: PETSC_INTERN PetscErrorCode EPSComputeVectors_Hermitian(EPS);
169: PETSC_INTERN PetscErrorCode EPSComputeVectors_Schur(EPS);
170: PETSC_INTERN PetscErrorCode EPSComputeVectors_Indefinite(EPS);
171: PETSC_INTERN PetscErrorCode EPSComputeVectors_Slice(EPS);
172: PETSC_INTERN PetscErrorCode EPSComputeResidualNorm_Private(EPS,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
173: PETSC_INTERN PetscErrorCode EPSComputeRitzVector(EPS,PetscScalar*,PetscScalar*,BV,Vec,Vec);
174: PETSC_INTERN PetscErrorCode EPSGetStartVector(EPS,PetscInt,PetscBool*);

176: /* Private functions of the solver implementations */

178: PETSC_INTERN PetscErrorCode EPSBasicArnoldi(EPS,PetscBool,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
179: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
180: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi1(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
181: PETSC_INTERN PetscErrorCode EPSKrylovConvergence(EPS,PetscBool,PetscInt,PetscInt,PetscReal,PetscReal,PetscInt*);
182: PETSC_INTERN PetscErrorCode EPSFullLanczos(EPS,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*);
183: PETSC_INTERN PetscErrorCode EPSPseudoLanczos(EPS,PetscReal*,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*,PetscBool*,PetscReal*,Vec);
184: PETSC_INTERN PetscErrorCode EPSBuildBalance_Krylov(EPS);

186: #endif