Actual source code: epsimpl.h

slepc-3.10.0 2018-09-18
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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: */

 11: #if !defined(_EPSIMPL)
 12: #define _EPSIMPL

 14: #include <slepceps.h>
 15: #include <slepc/private/slepcimpl.h>

 17: PETSC_EXTERN PetscBool EPSRegisterAllCalled;
 18: PETSC_EXTERN PetscErrorCode EPSRegisterAll(void);
 19: PETSC_EXTERN PetscLogEvent EPS_SetUp,EPS_Solve;

 21: typedef struct _EPSOps *EPSOps;

 23: struct _EPSOps {
 24:   PetscErrorCode (*solve)(EPS);
 25:   PetscErrorCode (*setup)(EPS);
 26:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,EPS);
 27:   PetscErrorCode (*publishoptions)(EPS);
 28:   PetscErrorCode (*destroy)(EPS);
 29:   PetscErrorCode (*reset)(EPS);
 30:   PetscErrorCode (*view)(EPS,PetscViewer);
 31:   PetscErrorCode (*backtransform)(EPS);
 32:   PetscErrorCode (*computevectors)(EPS);
 33:   PetscErrorCode (*setdefaultst)(EPS);
 34: };

 36: /*
 37:    Maximum number of monitors you can run with a single EPS
 38: */
 39: #define MAXEPSMONITORS 5

 41: /*
 42:    The solution process goes through several states
 43: */
 44: typedef enum { EPS_STATE_INITIAL,
 45:                EPS_STATE_SETUP,
 46:                EPS_STATE_SOLVED,
 47:                EPS_STATE_EIGENVECTORS } EPSStateType;

 49: /*
 50:    To classify the different solvers into categories
 51: */
 52: typedef enum { EPS_CATEGORY_KRYLOV,      /* Krylov solver: relies on STApply and STBackTransform (same as OTHER) */
 53:                EPS_CATEGORY_PRECOND,     /* Preconditioned solver: uses ST only to manage preconditioner */
 54:                EPS_CATEGORY_CONTOUR,     /* Contour integral: ST used to solve linear systems at integration points */
 55:                EPS_CATEGORY_OTHER } EPSSolverType;

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

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

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

117:   /* ----------------------- Status variables --------------------------*/
118:   EPSStateType   state;            /* initial -> setup -> solved -> eigenvectors */
119:   EPSSolverType  categ;            /* solver category */
120:   PetscInt       nconv;            /* number of converged eigenvalues */
121:   PetscInt       its;              /* number of iterations so far computed */
122:   PetscInt       n,nloc;           /* problem dimensions (global, local) */
123:   PetscReal      nrma,nrmb;        /* computed matrix norms */
124:   PetscBool      useds;            /* whether the solver uses the DS object or not */
125:   PetscBool      isgeneralized;
126:   PetscBool      ispositive;
127:   PetscBool      ishermitian;
128:   EPSConvergedReason reason;
129: };

131: /*
132:     Macros to test valid EPS arguments
133: */
134: #if !defined(PETSC_USE_DEBUG)

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

138: #else

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

145: #endif

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

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

167: /*
168:   EPS_Purify - purify the first k vectors in the V basis
169: */
170: PETSC_STATIC_INLINE PetscErrorCode EPS_Purify(EPS eps,PetscInt k)
171: {
173:   PetscInt       i;
174:   Vec            v,z;

177:   BVCreateVec(eps->V,&v);
178:   for (i=0;i<k;i++) {
179:     BVCopyVec(eps->V,i,v);
180:     BVGetColumn(eps->V,i,&z);
181:     STApply(eps->st,v,z);
182:     BVRestoreColumn(eps->V,i,&z);
183:   }
184:   VecDestroy(&v);
185:   return(0);
186: }

188: PETSC_INTERN PetscErrorCode EPSSetWhichEigenpairs_Default(EPS);
189: PETSC_INTERN PetscErrorCode EPSSetDimensions_Default(EPS,PetscInt,PetscInt*,PetscInt*);
190: PETSC_INTERN PetscErrorCode EPSBackTransform_Default(EPS);
191: PETSC_INTERN PetscErrorCode EPSComputeVectors(EPS);
192: PETSC_INTERN PetscErrorCode EPSComputeVectors_Hermitian(EPS);
193: PETSC_INTERN PetscErrorCode EPSComputeVectors_Schur(EPS);
194: PETSC_INTERN PetscErrorCode EPSComputeVectors_Indefinite(EPS);
195: PETSC_INTERN PetscErrorCode EPSComputeVectors_Slice(EPS);
196: PETSC_INTERN PetscErrorCode EPSComputeResidualNorm_Private(EPS,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
197: PETSC_INTERN PetscErrorCode EPSComputeRitzVector(EPS,PetscScalar*,PetscScalar*,BV,Vec,Vec);
198: PETSC_INTERN PetscErrorCode EPSGetStartVector(EPS,PetscInt,PetscBool*);

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

202: PETSC_INTERN PetscErrorCode EPSBasicArnoldi(EPS,PetscBool,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
203: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
204: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi1(EPS,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
205: PETSC_INTERN PetscErrorCode EPSKrylovConvergence(EPS,PetscBool,PetscInt,PetscInt,PetscReal,PetscReal,PetscInt*);
206: PETSC_INTERN PetscErrorCode EPSFullLanczos(EPS,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*);
207: PETSC_INTERN PetscErrorCode EPSPseudoLanczos(EPS,PetscReal*,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*,PetscBool*,PetscReal*,Vec);
208: PETSC_INTERN PetscErrorCode EPSBuildBalance_Krylov(EPS);
209: PETSC_INTERN PetscErrorCode EPSSetDefaultST(EPS);
210: PETSC_INTERN PetscErrorCode EPSSetDefaultST_Precond(EPS);
211: PETSC_INTERN PetscErrorCode EPSSetDefaultST_GMRES(EPS);

213: #endif