Actual source code: epsimpl.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(_EPSIMPL)
 23: #define _EPSIMPL

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

 28: PETSC_EXTERN PetscLogEvent EPS_SetUp,EPS_Solve;

 30: typedef struct _EPSOps *EPSOps;

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

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

 49: typedef enum { EPS_STATE_INITIAL,
 50:                EPS_STATE_SETUP,
 51:                EPS_STATE_SOLVED,
 52:                EPS_STATE_EIGENVECTORS } EPSStateType;

 54: /*
 55:    Defines the EPS data structure.
 56: */
 57: struct _p_EPS {
 58:   PETSCHEADER(struct _EPSOps);
 59:   /*------------------------- User parameters ---------------------------*/
 60:   PetscInt       max_it;           /* maximum number of iterations */
 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       nini;             /* number of initial vectors (negative means not copied yet) */
 65:   PetscInt       nds;              /* number of basis vectors of deflation space */
 66:   PetscScalar    target;           /* target value */
 67:   PetscReal      tol;              /* tolerance */
 68:   EPSConv        conv;             /* convergence test */
 69:   EPSWhich       which;            /* which part of the spectrum to be sought */
 70:   PetscReal      inta,intb;        /* interval [a,b] for spectrum slicing */
 71:   EPSProblemType problem_type;     /* which kind of problem to be solved */
 72:   EPSExtraction  extraction;       /* which kind of extraction to be applied */
 73:   EPSBalance     balance;          /* the balancing method */
 74:   PetscInt       balance_its;      /* number of iterations of the balancing method */
 75:   PetscReal      balance_cutoff;   /* cutoff value for balancing */
 76:   PetscBool      trueres;          /* whether the true residual norm must be computed */
 77:   PetscBool      trackall;         /* whether all the residuals must be computed */

 79:   /*-------------- User-provided functions and contexts -----------------*/
 80:   PetscErrorCode (*converged)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 81:   PetscErrorCode (*convergeddestroy)(void*);
 82:   PetscErrorCode (*arbitrary)(PetscScalar,PetscScalar,Vec,Vec,PetscScalar*,PetscScalar*,void*);
 83:   void           *convergedctx;
 84:   void           *arbitraryctx;
 85:   PetscErrorCode (*monitor[MAXEPSMONITORS])(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
 86:   PetscErrorCode (*monitordestroy[MAXEPSMONITORS])(void**);
 87:   void           *monitorcontext[MAXEPSMONITORS];
 88:   PetscInt       numbermonitors;

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

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

120: /*
121:     Macros to test valid EPS arguments
122: */
123: #if !defined(PETSC_USE_DEBUG)

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

127: #else

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

134: #endif

138: /*
139:   EPS_SetInnerProduct - set B matrix for inner product if appropriate.
140: */
141: PETSC_STATIC_INLINE PetscErrorCode EPS_SetInnerProduct(EPS eps)
142: {
144:   Mat            B;

147:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
148:   if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
149:     STGetBilinearForm(eps->st,&B);
150:     BVSetMatrix(eps->V,B,eps->ispositive?PETSC_FALSE:PETSC_TRUE);
151:     MatDestroy(&B);
152:   } else {
153:     BVSetMatrix(eps->V,NULL,PETSC_FALSE);
154:   }
155:   return(0);
156: }

158: PETSC_INTERN PetscErrorCode EPSSetWhichEigenpairs_Default(EPS);
159: PETSC_INTERN PetscErrorCode EPSSetDimensions_Default(EPS,PetscInt,PetscInt*,PetscInt*);
160: PETSC_INTERN PetscErrorCode EPSBackTransform_Default(EPS);
161: PETSC_INTERN PetscErrorCode EPSComputeVectors_Hermitian(EPS);
162: PETSC_INTERN PetscErrorCode EPSComputeVectors_Schur(EPS);
163: PETSC_INTERN PetscErrorCode EPSComputeVectors_Indefinite(EPS);
164: PETSC_INTERN PetscErrorCode EPSComputeResidualNorm_Private(EPS,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
165: PETSC_INTERN PetscErrorCode EPSComputeRelativeError_Private(EPS,PetscScalar,PetscScalar,Vec,Vec,PetscReal*);
166: PETSC_INTERN PetscErrorCode EPSComputeRitzVector(EPS,PetscScalar*,PetscScalar*,BV,Vec,Vec);
167: PETSC_INTERN PetscErrorCode EPSGetStartVector(EPS,PetscInt,PetscBool*);

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

171: PETSC_INTERN PetscErrorCode EPSBasicArnoldi(EPS,PetscBool,PetscScalar*,PetscInt,PetscInt,PetscInt*,PetscReal*,PetscBool*);
172: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi(EPS,PetscScalar*,PetscInt,Vec*,PetscInt,PetscInt*,Vec,PetscReal*,PetscBool*);
173: PETSC_INTERN PetscErrorCode EPSDelayedArnoldi1(EPS,PetscScalar*,PetscInt,Vec*,PetscInt,PetscInt*,Vec,PetscReal*,PetscBool*);
174: PETSC_INTERN PetscErrorCode EPSKrylovConvergence(EPS,PetscBool,PetscInt,PetscInt,PetscReal,PetscReal,PetscInt*);
175: PETSC_INTERN PetscErrorCode EPSFullLanczos(EPS,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*);
176: PETSC_INTERN PetscErrorCode EPSPseudoLanczos(EPS,PetscReal*,PetscReal*,PetscReal*,PetscInt,PetscInt*,PetscBool*,PetscReal*,Vec);
177: PETSC_INTERN PetscErrorCode EPSBuildBalance_Krylov(EPS);

179: #endif