Actual source code: pepimpl.h

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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(SLEPCPEPIMPL_H)
 12: #define SLEPCPEPIMPL_H

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

 17: SLEPC_EXTERN PetscBool PEPRegisterAllCalled;
 18: SLEPC_EXTERN PetscBool PEPMonitorRegisterAllCalled;
 19: SLEPC_EXTERN PetscErrorCode PEPRegisterAll(void);
 20: SLEPC_EXTERN PetscErrorCode PEPMonitorRegisterAll(void);
 21: SLEPC_EXTERN PetscLogEvent PEP_SetUp,PEP_Solve,PEP_Refine;

 23: typedef struct _PEPOps *PEPOps;

 25: struct _PEPOps {
 26:   PetscErrorCode (*solve)(PEP);
 27:   PetscErrorCode (*setup)(PEP);
 28:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PEP);
 29:   PetscErrorCode (*publishoptions)(PEP);
 30:   PetscErrorCode (*destroy)(PEP);
 31:   PetscErrorCode (*reset)(PEP);
 32:   PetscErrorCode (*view)(PEP,PetscViewer);
 33:   PetscErrorCode (*backtransform)(PEP);
 34:   PetscErrorCode (*computevectors)(PEP);
 35:   PetscErrorCode (*extractvectors)(PEP);
 36:   PetscErrorCode (*setdefaultst)(PEP);
 37: };

 39: /*
 40:      Maximum number of monitors you can run with a single PEP
 41: */
 42: #define MAXPEPMONITORS 5

 44: typedef enum { PEP_STATE_INITIAL,
 45:                PEP_STATE_SETUP,
 46:                PEP_STATE_SOLVED,
 47:                PEP_STATE_EIGENVECTORS } PEPStateType;

 49: /*
 50:    To check for unsupported features at PEPSetUp_XXX()
 51: */
 52: typedef enum { PEP_FEATURE_NONMONOMIAL=1,   /* non-monomial bases */
 53:                PEP_FEATURE_REGION=4,        /* nontrivial region for filtering */
 54:                PEP_FEATURE_EXTRACT=8,       /* eigenvector extraction */
 55:                PEP_FEATURE_CONVERGENCE=16,  /* convergence test selected by user */
 56:                PEP_FEATURE_STOPPING=32      /* stopping test */
 57:              } PEPFeatureType;

 59: /*
 60:    Defines the PEP data structure.
 61: */
 62: struct _p_PEP {
 63:   PETSCHEADER(struct _PEPOps);
 64:   /*------------------------- User parameters ---------------------------*/
 65:   PetscInt       max_it;           /* maximum number of iterations */
 66:   PetscInt       nev;              /* number of eigenvalues to compute */
 67:   PetscInt       ncv;              /* number of basis vectors */
 68:   PetscInt       mpd;              /* maximum dimension of projected problem */
 69:   PetscInt       nini;             /* number of initial vectors (negative means not copied yet) */
 70:   PetscScalar    target;           /* target value */
 71:   PetscReal      tol;              /* tolerance */
 72:   PEPConv        conv;             /* convergence test */
 73:   PEPStop        stop;             /* stopping test */
 74:   PEPWhich       which;            /* which part of the spectrum to be sought */
 75:   PetscReal      inta,intb;        /* interval [a,b] for spectrum slicing */
 76:   PEPBasis       basis;            /* polynomial basis used to represent the problem */
 77:   PEPProblemType problem_type;     /* which kind of problem to be solved */
 78:   PEPScale       scale;            /* scaling strategy to be used */
 79:   PetscReal      sfactor,dsfactor; /* scaling factors */
 80:   PetscInt       sits;             /* number of iterations of the scaling method */
 81:   PetscReal      slambda;          /* norm eigenvalue approximation for scaling */
 82:   PEPRefine      refine;           /* type of refinement to be applied after solve */
 83:   PetscInt       npart;            /* number of partitions of the communicator */
 84:   PetscReal      rtol;             /* tolerance for refinement */
 85:   PetscInt       rits;             /* number of iterations of the refinement method */
 86:   PEPRefineScheme scheme;          /* scheme for solving linear systems within refinement */
 87:   PEPExtract     extract;          /* type of extraction used */
 88:   PetscBool      trackall;         /* whether all the residuals must be computed */

 90:   /*-------------- User-provided functions and contexts -----------------*/
 91:   PetscErrorCode (*converged)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 92:   PetscErrorCode (*convergeduser)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 93:   PetscErrorCode (*convergeddestroy)(void*);
 94:   PetscErrorCode (*stopping)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
 95:   PetscErrorCode (*stoppinguser)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
 96:   PetscErrorCode (*stoppingdestroy)(void*);
 97:   void           *convergedctx;
 98:   void           *stoppingctx;
 99:   PetscErrorCode (*monitor[MAXPEPMONITORS])(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
100:   PetscErrorCode (*monitordestroy[MAXPEPMONITORS])(void**);
101:   void           *monitorcontext[MAXPEPMONITORS];
102:   PetscInt        numbermonitors;

104:   /*----------------- Child objects and working data -------------------*/
105:   ST             st;               /* spectral transformation object */
106:   DS             ds;               /* direct solver object */
107:   BV             V;                /* set of basis vectors and computed eigenvectors */
108:   RG             rg;               /* optional region for filtering */
109:   SlepcSC        sc;               /* sorting criterion data */
110:   Mat            *A;               /* coefficient matrices of the polynomial */
111:   PetscInt       nmat;             /* number of matrices */
112:   Vec            Dl,Dr;            /* diagonal matrices for balancing */
113:   Vec            *IS;              /* references to user-provided initial space */
114:   PetscScalar    *eigr,*eigi;      /* real and imaginary parts of eigenvalues */
115:   PetscReal      *errest;          /* error estimates */
116:   PetscInt       *perm;            /* permutation for eigenvalue ordering */
117:   PetscReal      *pbc;             /* coefficients defining the polynomial basis */
118:   PetscScalar    *solvematcoeffs;  /* coefficients to compute the matrix to be inverted */
119:   PetscInt       nwork;            /* number of work vectors */
120:   Vec            *work;            /* work vectors */
121:   KSP            refineksp;        /* ksp used in refinement */
122:   PetscSubcomm   refinesubc;       /* context for sub-communicators */
123:   void           *data;            /* placeholder for solver-specific stuff */

125:   /* ----------------------- Status variables --------------------------*/
126:   PEPStateType   state;            /* initial -> setup -> solved -> eigenvectors */
127:   PetscInt       nconv;            /* number of converged eigenvalues */
128:   PetscInt       its;              /* number of iterations so far computed */
129:   PetscInt       n,nloc;           /* problem dimensions (global, local) */
130:   PetscReal      *nrma;            /* computed matrix norms */
131:   PetscReal      nrml[2];          /* computed matrix norms for the linearization */
132:   PetscBool      sfactor_set;      /* flag to indicate the user gave sfactor */
133:   PetscBool      lineariz;         /* current solver is based on linearization */
134:   PEPConvergedReason reason;
135: };

137: /*
138:     Macros to test valid PEP arguments
139: */
140: #if !defined(PETSC_USE_DEBUG)

142: #define PEPCheckSolved(h,arg) do {(void)(h);} while (0)

144: #else

146: #define PEPCheckSolved(h,arg) \
147:   do { \
148:     if ((h)->state<PEP_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call PEPSolve() first: Parameter #%d",arg); \
149:   } while (0)

151: #endif

153: /*
154:     Macros to check settings at PEPSetUp()
155: */

157: /* PEPCheckHermitian: the problem is Hermitian or Hyperbolic */
158: #define PEPCheckHermitianCondition(pep,condition,msg) \
159:   do { \
160:     if (condition) { \
161:       if ((pep)->problem_type!=PEP_HERMITIAN && (pep)->problem_type!=PEP_HYPERBOLIC) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s can only be used for Hermitian (or hyperbolic) problems",((PetscObject)(pep))->type_name,(msg)); \
162:     } \
163:   } while (0)
164: #define PEPCheckHermitian(pep) PEPCheckHermitianCondition(pep,PETSC_TRUE,"")

166: /* PEPCheckQuadratic: the polynomial has degree 2 */
167: #define PEPCheckQuadraticCondition(pep,condition,msg) \
168:   do { \
169:     if (condition) { \
170:       if ((pep)->nmat!=3) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s is only available for quadratic problems",((PetscObject)(pep))->type_name,(msg)); \
171:     } \
172:   } while (0)
173: #define PEPCheckQuadratic(pep) PEPCheckQuadraticCondition(pep,PETSC_TRUE,"")

175: /* PEPCheckShiftSinvert: shift or shift-and-invert ST */
176: #define PEPCheckShiftSinvertCondition(pep,condition,msg) \
177:   do { \
178:     if (condition) { \
179:       PetscBool __flg; \
180:       PetscObjectTypeCompareAny((PetscObject)(pep)->st,&__flg,STSINVERT,STSHIFT,""); \
181:       if (!__flg) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s requires shift or shift-and-invert spectral transform",((PetscObject)(pep))->type_name,(msg)); \
182:     } \
183:   } while (0)
184: #define PEPCheckShiftSinvert(pep) PEPCheckShiftSinvertCondition(pep,PETSC_TRUE,"")

186: /* PEPCheckSinvertCayley: shift-and-invert or Cayley ST */
187: #define PEPCheckSinvertCayleyCondition(pep,condition,msg) \
188:   do { \
189:     if (condition) { \
190:       PetscBool __flg; \
191:       PetscObjectTypeCompareAny((PetscObject)(pep)->st,&__flg,STSINVERT,STCAYLEY,""); \
192:       if (!__flg) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s requires shift-and-invert or Cayley transform",((PetscObject)(pep))->type_name,(msg)); \
193:     } \
194:   } while (0)
195: #define PEPCheckSinvertCayley(pep) PEPCheckSinvertCayleyCondition(pep,PETSC_TRUE,"")

197: /* Check for unsupported features */
198: #define PEPCheckUnsupportedCondition(pep,mask,condition,msg) \
199:   do { \
200:     if (condition) { \
201:       if (((mask) & PEP_FEATURE_NONMONOMIAL) && (pep)->basis!=PEP_BASIS_MONOMIAL) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s is not implemented for non-monomial bases",((PetscObject)(pep))->type_name,(msg)); \
202:       if ((mask) & PEP_FEATURE_REGION) { \
203:         PetscBool      __istrivial; \
204:         PetscErrorCode __RGIsTrivial((pep)->rg,&__istrivial);CHKERRQ(__ierr); \
205:         if (!__istrivial) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s does not support region filtering",((PetscObject)(pep))->type_name,(msg)); \
206:       } \
207:       if (((mask) & PEP_FEATURE_EXTRACT) && (pep)->extract && (pep)->extract!=PEP_EXTRACT_NONE) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s does not support extraction variants",((PetscObject)(pep))->type_name,(msg)); \
208:       if (((mask) & PEP_FEATURE_CONVERGENCE) && (pep)->converged!=PEPConvergedRelative) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default convergence test",((PetscObject)(pep))->type_name,(msg)); \
209:       if (((mask) & PEP_FEATURE_STOPPING) && (pep)->stopping!=PEPStoppingBasic) SETERRQ2(PetscObjectComm((PetscObject)(pep)),PETSC_ERR_SUP,"The solver '%s'%s only supports the default stopping test",((PetscObject)(pep))->type_name,(msg)); \
210:     } \
211:   } while (0)
212: #define PEPCheckUnsupported(pep,mask) PEPCheckUnsupportedCondition(pep,mask,PETSC_TRUE,"")

214: /* Check for ignored features */
215: #define PEPCheckIgnoredCondition(pep,mask,condition,msg) \
216:   do { \
217:     PetscErrorCode __ierr; \
218:     if (condition) { \
219:       if (((mask) & PEP_FEATURE_NONMONOMIAL) && (pep)->basis!=PEP_BASIS_MONOMIAL) { __PetscInfo2((pep),"The solver '%s'%s ignores the basis settings\n",((PetscObject)(pep))->type_name,(msg)); } \
220:       if ((mask) & PEP_FEATURE_REGION) { \
221:         PetscBool __istrivial; \
222:         __RGIsTrivial((pep)->rg,&__istrivial);CHKERRQ(__ierr); \
223:         if (!__istrivial) { __PetscInfo2((pep),"The solver '%s'%s ignores the specified region\n",((PetscObject)(pep))->type_name,(msg)); } \
224:       } \
225:       if (((mask) & PEP_FEATURE_EXTRACT) && (pep)->extract && (pep)->extract!=PEP_EXTRACT_NONE) { __PetscInfo2((pep),"The solver '%s'%s ignores the extract settings\n",((PetscObject)(pep))->type_name,(msg)); } \
226:       if (((mask) & PEP_FEATURE_CONVERGENCE) && (pep)->converged!=PEPConvergedRelative) { __PetscInfo2((pep),"The solver '%s'%s ignores the convergence test settings\n",((PetscObject)(pep))->type_name,(msg)); } \
227:       if (((mask) & PEP_FEATURE_STOPPING) && (pep)->stopping!=PEPStoppingBasic) { __PetscInfo2((pep),"The solver '%s'%s ignores the stopping test settings\n",((PetscObject)(pep))->type_name,(msg)); } \
228:     } \
229:   } while (0)
230: #define PEPCheckIgnored(pep,mask) PEPCheckIgnoredCondition(pep,mask,PETSC_TRUE,"")


233: SLEPC_INTERN PetscErrorCode PEPSetWhichEigenpairs_Default(PEP);
234: SLEPC_INTERN PetscErrorCode PEPSetDimensions_Default(PEP,PetscInt,PetscInt*,PetscInt*);
235: SLEPC_INTERN PetscErrorCode PEPExtractVectors(PEP);
236: SLEPC_INTERN PetscErrorCode PEPBackTransform_Default(PEP);
237: SLEPC_INTERN PetscErrorCode PEPComputeVectors(PEP);
238: SLEPC_INTERN PetscErrorCode PEPComputeVectors_Default(PEP);
239: SLEPC_INTERN PetscErrorCode PEPComputeVectors_Indefinite(PEP);
240: SLEPC_INTERN PetscErrorCode PEPComputeResidualNorm_Private(PEP,PetscScalar,PetscScalar,Vec,Vec,Vec*,PetscReal*);
241: SLEPC_INTERN PetscErrorCode PEPKrylovConvergence(PEP,PetscBool,PetscInt,PetscInt,PetscReal,PetscInt*);
242: SLEPC_INTERN PetscErrorCode PEPComputeScaleFactor(PEP);
243: SLEPC_INTERN PetscErrorCode PEPBuildDiagonalScaling(PEP);
244: SLEPC_INTERN PetscErrorCode PEPBasisCoefficients(PEP,PetscReal*);
245: SLEPC_INTERN PetscErrorCode PEPEvaluateBasis(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
246: SLEPC_INTERN PetscErrorCode PEPEvaluateBasisDerivative(PEP,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
247: SLEPC_INTERN PetscErrorCode PEPEvaluateBasisMat(PEP,PetscInt,PetscScalar*,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt);
248: SLEPC_INTERN PetscErrorCode PEPNewtonRefinement_TOAR(PEP,PetscScalar,PetscInt*,PetscReal*,PetscInt,PetscScalar*,PetscInt);
249: SLEPC_INTERN PetscErrorCode PEPNewtonRefinementSimple(PEP,PetscInt*,PetscReal,PetscInt);
250: SLEPC_INTERN PetscErrorCode PEPSetDefaultST(PEP);
251: SLEPC_INTERN PetscErrorCode PEPSetDefaultST_Transform(PEP);

253: #endif