Actual source code: nepimpl.h

slepc-3.10.2 2019-02-11
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(_NEPIMPL)
 12: #define _NEPIMPL

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

 17: PETSC_EXTERN PetscBool NEPRegisterAllCalled;
 18: PETSC_EXTERN PetscErrorCode NEPRegisterAll(void);
 19: PETSC_EXTERN PetscLogEvent NEP_SetUp,NEP_Solve,NEP_Refine,NEP_FunctionEval,NEP_JacobianEval,NEP_DerivativesEval;

 21: typedef struct _NEPOps *NEPOps;

 23: struct _NEPOps {
 24:   PetscErrorCode (*solve)(NEP);
 25:   PetscErrorCode (*setup)(NEP);
 26:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,NEP);
 27:   PetscErrorCode (*publishoptions)(NEP);
 28:   PetscErrorCode (*destroy)(NEP);
 29:   PetscErrorCode (*reset)(NEP);
 30:   PetscErrorCode (*view)(NEP,PetscViewer);
 31:   PetscErrorCode (*computevectors)(NEP);
 32: };

 34: /*
 35:      Maximum number of monitors you can run with a single NEP
 36: */
 37: #define MAXNEPMONITORS 5

 39: typedef enum { NEP_STATE_INITIAL,
 40:                NEP_STATE_SETUP,
 41:                NEP_STATE_SOLVED,
 42:                NEP_STATE_EIGENVECTORS } NEPStateType;

 44: /*
 45:      How the problem function T(lambda) has been defined by the user
 46:      - Callback: one callback to build the function matrix, another one for the Jacobian
 47:      - Split: in split form sum_j(A_j*f_j(lambda))
 48:      - Derivatives: a single callback for all the derivatives (including the 0th derivative)
 49: */
 50: typedef enum { NEP_USER_INTERFACE_CALLBACK=1,
 51:                NEP_USER_INTERFACE_SPLIT,
 52:                NEP_USER_INTERFACE_DERIVATIVES } NEPUserInterface;

 54: /*
 55:    Defines the NEP data structure.
 56: */
 57: struct _p_NEP {
 58:   PETSCHEADER(struct _NEPOps);
 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:   PetscScalar    target;           /* target value */
 66:   PetscReal      tol;              /* tolerance */
 67:   NEPConv        conv;             /* convergence test */
 68:   NEPStop        stop;             /* stopping test */
 69:   NEPWhich       which;            /* which part of the spectrum to be sought */
 70:   NEPProblemType problem_type;     /* which kind of problem to be solved */
 71:   NEPRefine      refine;           /* type of refinement to be applied after solve */
 72:   PetscInt       npart;            /* number of partitions of the communicator */
 73:   PetscReal      rtol;             /* tolerance for refinement */
 74:   PetscInt       rits;             /* number of iterations of the refinement method */
 75:   NEPRefineScheme scheme;          /* scheme for solving linear systems within refinement */
 76:   PetscBool      trackall;         /* whether all the residuals must be computed */

 78:   /*-------------- User-provided functions and contexts -----------------*/
 79:   PetscErrorCode (*computefunction)(NEP,PetscScalar,Mat,Mat,void*);
 80:   PetscErrorCode (*computejacobian)(NEP,PetscScalar,Mat,void*);
 81:   void           *functionctx;
 82:   void           *jacobianctx;
 83:   PetscErrorCode (*computederivatives)(NEP,PetscScalar,PetscInt,Mat,void*);
 84:   void           *derivativesctx;
 85:   PetscErrorCode (*converged)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 86:   PetscErrorCode (*convergeduser)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
 87:   PetscErrorCode (*convergeddestroy)(void*);
 88:   PetscErrorCode (*stopping)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);
 89:   PetscErrorCode (*stoppinguser)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);
 90:   PetscErrorCode (*stoppingdestroy)(void*);
 91:   void           *convergedctx;
 92:   void           *stoppingctx;
 93:   PetscErrorCode (*monitor[MAXNEPMONITORS])(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
 94:   PetscErrorCode (*monitordestroy[MAXNEPMONITORS])(void**);
 95:   void           *monitorcontext[MAXNEPMONITORS];
 96:   PetscInt       numbermonitors;

 98:   /*----------------- Child objects and working data -------------------*/
 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:   Mat            function;         /* function matrix */
104:   Mat            function_pre;     /* function matrix (preconditioner) */
105:   Mat            jacobian;         /* Jacobian matrix */
106:   Mat            derivatives;      /* derivatives matrix */
107:   Mat            *A;               /* matrix coefficients of split form */
108:   FN             *f;               /* matrix functions of split form */
109:   PetscInt       nt;               /* number of terms in split form */
110:   MatStructure   mstr;             /* pattern of split matrices */
111:   Vec            *IS;              /* references to user-provided initial space */
112:   PetscScalar    *eigr,*eigi;      /* real and imaginary parts of eigenvalues */
113:   PetscReal      *errest;          /* error estimates */
114:   PetscInt       *perm;            /* permutation for eigenvalue ordering */
115:   PetscInt       nwork;            /* number of work vectors */
116:   Vec            *work;            /* work vectors */
117:   KSP            refineksp;        /* ksp used in refinement */
118:   PetscSubcomm   refinesubc;       /* context for sub-communicators */
119:   void           *data;            /* placeholder for solver-specific stuff */

121:   /* ----------------------- Status variables --------------------------*/
122:   NEPStateType   state;            /* initial -> setup -> solved -> eigenvectors */
123:   PetscInt       nconv;            /* number of converged eigenvalues */
124:   PetscInt       its;              /* number of iterations so far computed */
125:   PetscInt       n,nloc;           /* problem dimensions (global, local) */
126:   PetscReal      *nrma;            /* computed matrix norms */
127:   NEPUserInterface fui;            /* how the user has defined the nonlinear operator */
128:   NEPConvergedReason reason;
129: };

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

136: #define NEPCheckProblem(h,arg) do {} while (0)
137: #define NEPCheckCallback(h,arg) do {} while (0)
138: #define NEPCheckSplit(h,arg) do {} while (0)
139: #define NEPCheckDerivatives(h,arg) do {} while (0)
140: #define NEPCheckSolved(h,arg) do {} while (0)

142: #else

144: #define NEPCheckProblem(h,arg) \
145:   do { \
146:     if (!(h->fui)) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"The nonlinear eigenproblem has not been specified yet. Parameter #%d",arg); \
147:   } while (0)

149: #define NEPCheckCallback(h,arg) \
150:   do { \
151:     if (h->fui!=NEP_USER_INTERFACE_CALLBACK) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"This operation requires the nonlinear eigenproblem specified with callbacks. Parameter #%d",arg); \
152:   } while (0)

154: #define NEPCheckSplit(h,arg) \
155:   do { \
156:     if (h->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"This operation requires the nonlinear eigenproblem in split form. Parameter #%d",arg); \
157:   } while (0)

159: #define NEPCheckDerivatives(h,arg) \
160:   do { \
161:     if (h->fui!=NEP_USER_INTERFACE_DERIVATIVES) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"This operation requires the nonlinear eigenproblem specified with derivatives callback. Parameter #%d",arg); \
162:   } while (0)

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

169: #endif

171: PETSC_INTERN PetscErrorCode NEPSetDimensions_Default(NEP,PetscInt,PetscInt*,PetscInt*);
172: PETSC_INTERN PetscErrorCode NEPComputeVectors(NEP);
173: PETSC_INTERN PetscErrorCode NEPReset_Problem(NEP);
174: PETSC_INTERN PetscErrorCode NEPGetDefaultShift(NEP,PetscScalar*);
175: PETSC_INTERN PetscErrorCode NEPComputeVectors_Schur(NEP);
176: PETSC_INTERN PetscErrorCode NEPComputeResidualNorm_Private(NEP,PetscScalar,Vec,Vec*,PetscReal*);
177: PETSC_INTERN PetscErrorCode NEPNewtonRefinementSimple(NEP,PetscInt*,PetscReal,PetscInt);

179: #endif