Actual source code: nepimpl.h

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

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

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

 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: */
 49: typedef enum { NEP_USER_INTERFACE_CALLBACK=1,
 50:                NEP_USER_INTERFACE_SPLIT } NEPUserInterface;

 52: /*
 53:    Defines the NEP data structure.
 54: */
 55: struct _p_NEP {
 56:   PETSCHEADER(struct _NEPOps);
 57:   /*------------------------- User parameters ---------------------------*/
 58:   PetscInt       max_it;           /* maximum number of iterations */
 59:   PetscInt       nev;              /* number of eigenvalues to compute */
 60:   PetscInt       ncv;              /* number of basis vectors */
 61:   PetscInt       mpd;              /* maximum dimension of projected problem */
 62:   PetscInt       nini;             /* number of initial vectors (negative means not copied yet) */
 63:   PetscScalar    target;           /* target value */
 64:   PetscReal      tol;              /* tolerance */
 65:   NEPConv        conv;             /* convergence test */
 66:   NEPStop        stop;             /* stopping test */
 67:   NEPWhich       which;            /* which part of the spectrum to be sought */
 68:   NEPProblemType problem_type;     /* which kind of problem to be solved */
 69:   NEPRefine      refine;           /* type of refinement to be applied after solve */
 70:   PetscInt       npart;            /* number of partitions of the communicator */
 71:   PetscReal      rtol;             /* tolerance for refinement */
 72:   PetscInt       rits;             /* number of iterations of the refinement method */
 73:   NEPRefineScheme scheme;          /* scheme for solving linear systems within refinement */
 74:   PetscBool      trackall;         /* whether all the residuals must be computed */
 75:   PetscBool      twosided;         /* whether to compute left eigenvectors (two-sided solver) */

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

 95:   /*----------------- Child objects and working data -------------------*/
 96:   DS             ds;               /* direct solver object */
 97:   BV             V;                /* set of basis vectors and computed eigenvectors */
 98:   BV             W;                /* left basis vectors (if left eigenvectors requested) */
 99:   RG             rg;               /* optional region for filtering */
100:   SlepcSC        sc;               /* sorting criterion data */
101:   Mat            function;         /* function matrix */
102:   Mat            function_pre;     /* function matrix (preconditioner) */
103:   Mat            jacobian;         /* Jacobian matrix */
104:   Mat            *A;               /* matrix coefficients of split form */
105:   FN             *f;               /* matrix functions of split form */
106:   PetscInt       nt;               /* number of terms in split form */
107:   MatStructure   mstr;             /* pattern of split matrices */
108:   Vec            *IS;              /* references to user-provided initial space */
109:   PetscScalar    *eigr,*eigi;      /* real and imaginary parts of eigenvalues */
110:   PetscReal      *errest;          /* error estimates */
111:   PetscInt       *perm;            /* permutation for eigenvalue ordering */
112:   PetscInt       nwork;            /* number of work vectors */
113:   Vec            *work;            /* work vectors */
114:   KSP            refineksp;        /* ksp used in refinement */
115:   PetscSubcomm   refinesubc;       /* context for sub-communicators */
116:   void           *data;            /* placeholder for solver-specific stuff */

118:   /* ----------------------- Status variables --------------------------*/
119:   NEPStateType   state;            /* initial -> setup -> solved -> eigenvectors */
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;            /* computed matrix norms */
124:   NEPUserInterface fui;            /* how the user has defined the nonlinear operator */
125:   PetscBool      useds;            /* whether the solver uses the DS object or not */
126:   PetscBool      hasts;            /* whether the solver has two-sided variant */
127:   Mat            resolvent;        /* shell matrix to be used in NEPApplyResolvent */
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 NEPCheckSolved(h,arg) \
160:   do { \
161:     if ((h)->state<NEP_STATE_SOLVED) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"Must call NEPSolve() first: Parameter #%d",arg); \
162:   } while (0)

164: #endif

166: SLEPC_INTERN PetscErrorCode NEPSetDimensions_Default(NEP,PetscInt,PetscInt*,PetscInt*);
167: SLEPC_INTERN PetscErrorCode NEPComputeVectors(NEP);
168: SLEPC_INTERN PetscErrorCode NEPReset_Problem(NEP);
169: SLEPC_INTERN PetscErrorCode NEPGetDefaultShift(NEP,PetscScalar*);
170: SLEPC_INTERN PetscErrorCode NEPComputeVectors_Schur(NEP);
171: SLEPC_INTERN PetscErrorCode NEPComputeResidualNorm_Private(NEP,PetscBool,PetscScalar,Vec,Vec*,PetscReal*);
172: SLEPC_INTERN PetscErrorCode NEPNewtonRefinementSimple(NEP,PetscInt*,PetscReal,PetscInt);

174: #endif