Actual source code: slepcnep.h

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: */
 10: /*
 11:    User interface for SLEPc's nonlinear eigenvalue solvers
 12: */

 14: #pragma once

 16: #include <slepceps.h>
 17: #include <slepcpep.h>
 18: #include <slepcfn.h>

 20: /* SUBMANSEC = NEP */

 22: SLEPC_EXTERN PetscErrorCode NEPInitializePackage(void);
 23: SLEPC_EXTERN PetscErrorCode NEPFinalizePackage(void);

 25: /*S
 26:      NEP - Abstract SLEPc object that manages all solvers for
 27:      nonlinear eigenvalue problems.

 29:    Level: beginner

 31: .seealso:  NEPCreate()
 32: S*/
 33: typedef struct _p_NEP* NEP;

 35: /*J
 36:     NEPType - String with the name of a nonlinear eigensolver

 38:    Level: beginner

 40: .seealso: NEPSetType(), NEP
 41: J*/
 42: typedef const char *NEPType;
 43: #define NEPRII       "rii"
 44: #define NEPSLP       "slp"
 45: #define NEPNARNOLDI  "narnoldi"
 46: #define NEPCISS      "ciss"
 47: #define NEPINTERPOL  "interpol"
 48: #define NEPNLEIGS    "nleigs"

 50: /* Logging support */
 51: SLEPC_EXTERN PetscClassId NEP_CLASSID;

 53: /*E
 54:     NEPProblemType - Determines the type of the nonlinear eigenproblem

 56:     Level: intermediate

 58: .seealso: NEPSetProblemType(), NEPGetProblemType()
 59: E*/
 60: typedef enum { NEP_GENERAL  = 1,
 61:                NEP_RATIONAL = 2     /* NEP defined in split form with all f_i rational */
 62:              } NEPProblemType;

 64: /*E
 65:     NEPWhich - Determines which part of the spectrum is requested

 67:     Level: intermediate

 69: .seealso: NEPSetWhichEigenpairs(), NEPGetWhichEigenpairs()
 70: E*/
 71: typedef enum { NEP_LARGEST_MAGNITUDE  = 1,
 72:                NEP_SMALLEST_MAGNITUDE = 2,
 73:                NEP_LARGEST_REAL       = 3,
 74:                NEP_SMALLEST_REAL      = 4,
 75:                NEP_LARGEST_IMAGINARY  = 5,
 76:                NEP_SMALLEST_IMAGINARY = 6,
 77:                NEP_TARGET_MAGNITUDE   = 7,
 78:                NEP_TARGET_REAL        = 8,
 79:                NEP_TARGET_IMAGINARY   = 9,
 80:                NEP_ALL                = 10,
 81:                NEP_WHICH_USER         = 11 } NEPWhich;

 83: /*E
 84:     NEPErrorType - The error type used to assess accuracy of computed solutions

 86:     Level: intermediate

 88: .seealso: NEPComputeError()
 89: E*/
 90: typedef enum { NEP_ERROR_ABSOLUTE,
 91:                NEP_ERROR_RELATIVE,
 92:                NEP_ERROR_BACKWARD } NEPErrorType;
 93: SLEPC_EXTERN const char *NEPErrorTypes[];

 95: /*E
 96:     NEPRefine - The refinement type

 98:     Level: intermediate

100: .seealso: NEPSetRefine()
101: E*/
102: typedef enum { NEP_REFINE_NONE,
103:                NEP_REFINE_SIMPLE,
104:                NEP_REFINE_MULTIPLE } NEPRefine;
105: SLEPC_EXTERN const char *NEPRefineTypes[];

107: /*E
108:     NEPRefineScheme - The scheme used for solving linear systems during iterative refinement

110:     Level: intermediate

112: .seealso: NEPSetRefine()
113: E*/
114: typedef enum { NEP_REFINE_SCHEME_SCHUR    = 1,
115:                NEP_REFINE_SCHEME_MBE      = 2,
116:                NEP_REFINE_SCHEME_EXPLICIT = 3 } NEPRefineScheme;
117: SLEPC_EXTERN const char *NEPRefineSchemes[];

119: /*E
120:     NEPConv - Determines the convergence test

122:     Level: intermediate

124: .seealso: NEPSetConvergenceTest(), NEPSetConvergenceTestFunction()
125: E*/
126: typedef enum { NEP_CONV_ABS,
127:                NEP_CONV_REL,
128:                NEP_CONV_NORM,
129:                NEP_CONV_USER } NEPConv;

131: /*E
132:     NEPStop - Determines the stopping test

134:     Level: advanced

136: .seealso: NEPSetStoppingTest(), NEPSetStoppingTestFunction()
137: E*/
138: typedef enum { NEP_STOP_BASIC,
139:                NEP_STOP_USER } NEPStop;

141: /*E
142:     NEPConvergedReason - Reason a nonlinear eigensolver was said to
143:          have converged or diverged

145:     Level: intermediate

147: .seealso: NEPSolve(), NEPGetConvergedReason(), NEPSetTolerances()
148: E*/
149: typedef enum {/* converged */
150:               NEP_CONVERGED_TOL                =  1,
151:               NEP_CONVERGED_USER               =  2,
152:               /* diverged */
153:               NEP_DIVERGED_ITS                 = -1,
154:               NEP_DIVERGED_BREAKDOWN           = -2,
155:                     /* unused                  = -3 */
156:               NEP_DIVERGED_LINEAR_SOLVE        = -4,
157:               NEP_DIVERGED_SUBSPACE_EXHAUSTED  = -5,
158:               NEP_CONVERGED_ITERATING          =  0} NEPConvergedReason;
159: SLEPC_EXTERN const char *const*NEPConvergedReasons;

161: SLEPC_EXTERN PetscErrorCode NEPCreate(MPI_Comm,NEP*);
162: SLEPC_EXTERN PetscErrorCode NEPDestroy(NEP*);
163: SLEPC_EXTERN PetscErrorCode NEPReset(NEP);
164: SLEPC_EXTERN PetscErrorCode NEPSetType(NEP,NEPType);
165: SLEPC_EXTERN PetscErrorCode NEPGetType(NEP,NEPType*);
166: SLEPC_EXTERN PetscErrorCode NEPSetProblemType(NEP,NEPProblemType);
167: SLEPC_EXTERN PetscErrorCode NEPGetProblemType(NEP,NEPProblemType*);
168: SLEPC_EXTERN PetscErrorCode NEPSetTarget(NEP,PetscScalar);
169: SLEPC_EXTERN PetscErrorCode NEPGetTarget(NEP,PetscScalar*);
170: SLEPC_EXTERN PetscErrorCode NEPSetFromOptions(NEP);
171: SLEPC_EXTERN PetscErrorCode NEPSetDSType(NEP);
172: SLEPC_EXTERN PetscErrorCode NEPSetUp(NEP);
173: SLEPC_EXTERN PetscErrorCode NEPSolve(NEP);
174: SLEPC_EXTERN PetscErrorCode NEPView(NEP,PetscViewer);
175: SLEPC_EXTERN PetscErrorCode NEPViewFromOptions(NEP,PetscObject,const char[]);
176: SLEPC_EXTERN PetscErrorCode NEPErrorView(NEP,NEPErrorType,PetscViewer);
177: SLEPC_EXTERN PetscErrorCode NEPErrorViewFromOptions(NEP);
178: SLEPC_EXTERN PetscErrorCode NEPConvergedReasonView(NEP,PetscViewer);
179: SLEPC_EXTERN PetscErrorCode NEPConvergedReasonViewFromOptions(NEP);
180: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "NEPConvergedReasonView()", ) static inline PetscErrorCode NEPReasonView(NEP nep,PetscViewer v) {return NEPConvergedReasonView(nep,v);}
181: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "NEPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode NEPReasonViewFromOptions(NEP nep) {return NEPConvergedReasonViewFromOptions(nep);}
182: SLEPC_EXTERN PetscErrorCode NEPValuesView(NEP,PetscViewer);
183: SLEPC_EXTERN PetscErrorCode NEPValuesViewFromOptions(NEP);
184: SLEPC_EXTERN PetscErrorCode NEPVectorsView(NEP,PetscViewer);
185: SLEPC_EXTERN PetscErrorCode NEPVectorsViewFromOptions(NEP);

187: /*S
188:   NEPFunctionFn - A prototype of a NEP function evaluation function that would be passed to NEPSetFunction()

190:   Calling Sequence:
191: +   nep    - eigensolver context obtained from NEPCreate()
192: .   lambda - the scalar argument where T(.) must be evaluated
193: .   T      - matrix that will contain T(lambda)
194: .   P      - [optional] different matrix to build the preconditioner
195: -   ctx    - [optional] user-defined context for private data for the
196:              function evaluation routine (may be NULL)

198:   Level: beginner

200: .seealso: NEPSetFunction()
201: S*/
202: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPFunctionFn(NEP nep,PetscScalar lambda,Mat T,Mat P,void *ctx);

204: /*S
205:   NEPJacobianFn - A prototype of a NEP Jacobian evaluation function that would be passed to NEPSetJacobian()

207:   Calling Sequence:
208: +   nep    - eigensolver context obtained from NEPCreate()
209: .   lambda - the scalar argument where T'(.) must be evaluated
210: .   J      - matrix that will contain T'(lambda)
211: -   ctx    - [optional] user-defined context for private data for the
212:              Jacobian evaluation routine (may be NULL)

214:   Level: beginner

216: .seealso: NEPSetJacobian()
217: S*/
218: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPJacobianFn(NEP nep,PetscScalar lambda,Mat J,void *ctx);

220: SLEPC_EXTERN PetscErrorCode NEPSetFunction(NEP,Mat,Mat,NEPFunctionFn*,void*);
221: SLEPC_EXTERN PetscErrorCode NEPGetFunction(NEP,Mat*,Mat*,NEPFunctionFn**,void**);
222: SLEPC_EXTERN PetscErrorCode NEPSetJacobian(NEP,Mat,NEPJacobianFn*,void*);
223: SLEPC_EXTERN PetscErrorCode NEPGetJacobian(NEP,Mat*,NEPJacobianFn**,void**);
224: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "NEPSetFunction() and NEPSetJacobian()", ) static inline PetscErrorCode NEPSetDerivatives(NEP nep,Mat A,PetscErrorCode (*fun)(NEP,PetscScalar,PetscInt,Mat,void*),void *ctx) {(void)A;(void)fun;(void)ctx;SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented in this version");}
225: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "NEPGetFunction() and NEPGetJacobian()", ) static inline PetscErrorCode NEPGetDerivatives(NEP nep,Mat *A,PetscErrorCode (**fun)(NEP,PetscScalar,PetscInt,Mat,void*),void **ctx) {(void)A;(void)fun;(void)ctx;SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented in this version");}
226: SLEPC_EXTERN PetscErrorCode NEPSetSplitOperator(NEP,PetscInt,Mat[],FN[],MatStructure);
227: SLEPC_EXTERN PetscErrorCode NEPGetSplitOperatorTerm(NEP,PetscInt,Mat*,FN*);
228: SLEPC_EXTERN PetscErrorCode NEPGetSplitOperatorInfo(NEP,PetscInt*,MatStructure*);
229: SLEPC_EXTERN PetscErrorCode NEPSetSplitPreconditioner(NEP,PetscInt,Mat[],MatStructure);
230: SLEPC_EXTERN PetscErrorCode NEPGetSplitPreconditionerTerm(NEP,PetscInt,Mat*);
231: SLEPC_EXTERN PetscErrorCode NEPGetSplitPreconditionerInfo(NEP,PetscInt*,MatStructure*);

233: SLEPC_EXTERN PetscErrorCode NEPSetBV(NEP,BV);
234: SLEPC_EXTERN PetscErrorCode NEPGetBV(NEP,BV*);
235: SLEPC_EXTERN PetscErrorCode NEPSetRG(NEP,RG);
236: SLEPC_EXTERN PetscErrorCode NEPGetRG(NEP,RG*);
237: SLEPC_EXTERN PetscErrorCode NEPSetDS(NEP,DS);
238: SLEPC_EXTERN PetscErrorCode NEPGetDS(NEP,DS*);
239: SLEPC_EXTERN PetscErrorCode NEPRefineGetKSP(NEP,KSP*);
240: SLEPC_EXTERN PetscErrorCode NEPSetTolerances(NEP,PetscReal,PetscInt);
241: SLEPC_EXTERN PetscErrorCode NEPGetTolerances(NEP,PetscReal*,PetscInt*);
242: SLEPC_EXTERN PetscErrorCode NEPSetDimensions(NEP,PetscInt,PetscInt,PetscInt);
243: SLEPC_EXTERN PetscErrorCode NEPGetDimensions(NEP,PetscInt*,PetscInt*,PetscInt*);
244: SLEPC_EXTERN PetscErrorCode NEPSetRefine(NEP,NEPRefine,PetscInt,PetscReal,PetscInt,NEPRefineScheme);
245: SLEPC_EXTERN PetscErrorCode NEPGetRefine(NEP,NEPRefine*,PetscInt*,PetscReal*,PetscInt*,NEPRefineScheme*);

247: SLEPC_EXTERN PetscErrorCode NEPGetConverged(NEP,PetscInt*);
248: SLEPC_EXTERN PetscErrorCode NEPGetEigenpair(NEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
249: SLEPC_EXTERN PetscErrorCode NEPGetLeftEigenvector(NEP,PetscInt,Vec,Vec);

251: SLEPC_EXTERN PetscErrorCode NEPComputeError(NEP,PetscInt,NEPErrorType,PetscReal*);
252: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "NEPComputeError()", ) static inline PetscErrorCode NEPComputeRelativeError(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_RELATIVE,r);}
253: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "NEPComputeError() with NEP_ERROR_ABSOLUTE", ) static inline PetscErrorCode NEPComputeResidualNorm(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_ABSOLUTE,r);}
254: SLEPC_EXTERN PetscErrorCode NEPGetErrorEstimate(NEP,PetscInt,PetscReal*);

256: SLEPC_EXTERN PetscErrorCode NEPComputeFunction(NEP,PetscScalar,Mat,Mat);
257: SLEPC_EXTERN PetscErrorCode NEPComputeJacobian(NEP,PetscScalar,Mat);
258: SLEPC_EXTERN PetscErrorCode NEPApplyFunction(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
259: SLEPC_EXTERN PetscErrorCode NEPApplyAdjoint(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
260: SLEPC_EXTERN PetscErrorCode NEPApplyJacobian(NEP,PetscScalar,Vec,Vec,Vec,Mat);
261: SLEPC_EXTERN PetscErrorCode NEPProjectOperator(NEP,PetscInt,PetscInt);
262: SLEPC_EXTERN PetscErrorCode NEPGetIterationNumber(NEP,PetscInt*);

264: SLEPC_EXTERN PetscErrorCode NEPSetInitialSpace(NEP,PetscInt,Vec[]);
265: SLEPC_EXTERN PetscErrorCode NEPSetWhichEigenpairs(NEP,NEPWhich);
266: SLEPC_EXTERN PetscErrorCode NEPGetWhichEigenpairs(NEP,NEPWhich*);
267: SLEPC_EXTERN PetscErrorCode NEPSetTwoSided(NEP,PetscBool);
268: SLEPC_EXTERN PetscErrorCode NEPGetTwoSided(NEP,PetscBool*);

270: SLEPC_EXTERN PetscErrorCode NEPApplyResolvent(NEP,RG,PetscScalar,Vec,Vec);

272: SLEPC_EXTERN PetscErrorCode NEPSetTrackAll(NEP,PetscBool);
273: SLEPC_EXTERN PetscErrorCode NEPGetTrackAll(NEP,PetscBool*);

275: SLEPC_EXTERN PetscErrorCode NEPGetConvergedReason(NEP,NEPConvergedReason*);

277: /*S
278:   NEPMonitorFn - A function prototype for functions provided to NEPMonitorSet()

280:   Calling Sequence:
281: +   nep    - eigensolver context obtained from NEPCreate()
282: .   its    - iteration number
283: .   nconv  - number of converged eigenpairs
284: .   eigr   - real part of the eigenvalues
285: .   eigi   - imaginary part of the eigenvalues
286: .   errest - relative error estimates for each eigenpair
287: .   nest   - number of error estimates
288: -   ctx    - optional monitoring context, as provided with NEPMonitorSet()

290:   Level: beginner

292: .seealso: NEPMonitorSet()
293: S*/
294: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPMonitorFn(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx);

296: /*S
297:   NEPMonitorRegisterFn - A function prototype for functions provided to NEPMonitorRegister()

299:   Calling Sequence:
300: +   nep    - eigensolver context obtained from NEPCreate()
301: .   its    - iteration number
302: .   nconv  - number of converged eigenpairs
303: .   eigr   - real part of the eigenvalues
304: .   eigi   - imaginary part of the eigenvalues
305: .   errest - relative error estimates for each eigenpair
306: .   nest   - number of error estimates
307: -   ctx    - PetscViewerAndFormat object

309:   Level: beginner

311:   Note:
312:   This is an NEPMonitorFn specialized for a context of PetscViewerAndFormat.

314: .seealso: NEPMonitorSet(), NEPMonitorRegister(), NEPMonitorFn, NEPMonitorRegisterCreateFn, NEPMonitorRegisterDestroyFn
315: S*/
316: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPMonitorRegisterFn(NEP nep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *ctx);

318: /*S
319:   NEPMonitorRegisterCreateFn - A function prototype for functions that do the creation when provided to NEPMonitorRegister()

321:   Calling Sequence:
322: +   viewer - the viewer to be used with the NEPMonitorRegisterFn
323: .   format - the format of the viewer
324: .   ctx    - a context for the monitor
325: -   result - a PetscViewerAndFormat object

327:   Level: beginner

329: .seealso: NEPMonitorRegisterFn, NEPMonitorSet(), NEPMonitorRegister(), NEPMonitorFn, NEPMonitorRegisterDestroyFn
330: S*/
331: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPMonitorRegisterCreateFn(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **result);

333: /*S
334:   NEPMonitorRegisterDestroyFn - A function prototype for functions that do the after use destruction when provided to NEPMonitorRegister()

336:   Calling Sequence:
337: .   vf - a PetscViewerAndFormat object to be destroyed, including any context

339:   Level: beginner

341: .seealso: NEPMonitorRegisterFn, NEPMonitorSet(), NEPMonitorRegister(), NEPMonitorFn, NEPMonitorRegisterCreateFn
342: S*/
343: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPMonitorRegisterDestroyFn(PetscViewerAndFormat **result);

345: SLEPC_EXTERN PetscErrorCode NEPMonitor(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
346: SLEPC_EXTERN PetscErrorCode NEPMonitorSet(NEP,NEPMonitorFn,void*,PetscCtxDestroyFn*);
347: SLEPC_EXTERN PetscErrorCode NEPMonitorCancel(NEP);
348: SLEPC_EXTERN PetscErrorCode NEPGetMonitorContext(NEP,void*);

350: SLEPC_EXTERN PetscErrorCode NEPMonitorSetFromOptions(NEP,const char[],const char[],void*,PetscBool);
351: SLEPC_EXTERN NEPMonitorRegisterFn        NEPMonitorFirst;
352: SLEPC_EXTERN NEPMonitorRegisterFn        NEPMonitorFirstDrawLG;
353: SLEPC_EXTERN NEPMonitorRegisterCreateFn  NEPMonitorFirstDrawLGCreate;
354: SLEPC_EXTERN NEPMonitorRegisterFn        NEPMonitorAll;
355: SLEPC_EXTERN NEPMonitorRegisterFn        NEPMonitorAllDrawLG;
356: SLEPC_EXTERN NEPMonitorRegisterCreateFn  NEPMonitorAllDrawLGCreate;
357: SLEPC_EXTERN NEPMonitorRegisterFn        NEPMonitorConverged;
358: SLEPC_EXTERN NEPMonitorRegisterCreateFn  NEPMonitorConvergedCreate;
359: SLEPC_EXTERN NEPMonitorRegisterFn        NEPMonitorConvergedDrawLG;
360: SLEPC_EXTERN NEPMonitorRegisterCreateFn  NEPMonitorConvergedDrawLGCreate;
361: SLEPC_EXTERN NEPMonitorRegisterDestroyFn NEPMonitorConvergedDestroy;

363: SLEPC_EXTERN PetscErrorCode NEPSetOptionsPrefix(NEP,const char*);
364: SLEPC_EXTERN PetscErrorCode NEPAppendOptionsPrefix(NEP,const char*);
365: SLEPC_EXTERN PetscErrorCode NEPGetOptionsPrefix(NEP,const char*[]);

367: SLEPC_EXTERN PetscFunctionList NEPList;
368: SLEPC_EXTERN PetscFunctionList NEPMonitorList;
369: SLEPC_EXTERN PetscFunctionList NEPMonitorCreateList;
370: SLEPC_EXTERN PetscFunctionList NEPMonitorDestroyList;
371: SLEPC_EXTERN PetscErrorCode NEPRegister(const char[],PetscErrorCode(*)(NEP));
372: SLEPC_EXTERN PetscErrorCode NEPMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,NEPMonitorRegisterFn*,NEPMonitorRegisterCreateFn*,NEPMonitorRegisterDestroyFn*);

374: SLEPC_EXTERN PetscErrorCode NEPSetWorkVecs(NEP,PetscInt);
375: SLEPC_EXTERN PetscErrorCode NEPAllocateSolution(NEP,PetscInt);

377: /*S
378:   NEPConvergenceTestFn - A prototype of a NEP convergence test function that would be passed to NEPSetConvergenceTestFunction()

380:   Calling Sequence:
381: +   nep    - eigensolver context obtained from NEPCreate()
382: .   eigr   - real part of the eigenvalue
383: .   eigi   - imaginary part of the eigenvalue
384: .   res    - residual norm associated to the eigenpair
385: .   errest - [output] computed error estimate
386: -   ctx    - [optional] user-defined context for private data for the
387:              convergence test routine (may be NULL)

389:   Level: advanced

391: .seealso: NEPSetConvergenceTestFunction()
392: S*/
393: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPConvergenceTestFn(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx);

395: SLEPC_EXTERN PetscErrorCode NEPSetConvergenceTest(NEP,NEPConv);
396: SLEPC_EXTERN PetscErrorCode NEPGetConvergenceTest(NEP,NEPConv*);
397: SLEPC_EXTERN NEPConvergenceTestFn NEPConvergedAbsolute;
398: SLEPC_EXTERN NEPConvergenceTestFn NEPConvergedRelative;
399: SLEPC_EXTERN NEPConvergenceTestFn NEPConvergedNorm;
400: SLEPC_EXTERN PetscErrorCode NEPSetConvergenceTestFunction(NEP,NEPConvergenceTestFn*,void*,PetscCtxDestroyFn*);

402: /*S
403:   NEPStoppingTestFn - A prototype of a NEP stopping test function that would be passed to NEPSetStoppingTestFunction()

405:   Calling Sequence:
406: +   nep    - eigensolver context obtained from NEPCreate()
407: .   its    - current number of iterations
408: .   max_it - maximum number of iterations
409: .   nconv  - number of currently converged eigenpairs
410: .   nev    - number of requested eigenpairs
411: .   reason - [output] result of the stopping test
412: -   ctx    - [optional] user-defined context for private data for the
413:              stopping test routine (may be NULL)

415:   Level: advanced

417: .seealso: NEPSetStoppingTestFunction()
418: S*/
419: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPStoppingTestFn(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx);

421: SLEPC_EXTERN PetscErrorCode NEPSetStoppingTest(NEP,NEPStop);
422: SLEPC_EXTERN PetscErrorCode NEPGetStoppingTest(NEP,NEPStop*);
423: SLEPC_EXTERN NEPStoppingTestFn NEPStoppingBasic;
424: SLEPC_EXTERN PetscErrorCode NEPSetStoppingTestFunction(NEP,NEPStoppingTestFn*,void*,PetscCtxDestroyFn*);

426: SLEPC_EXTERN PetscErrorCode NEPSetEigenvalueComparison(NEP,SlepcEigenvalueComparisonFn*,void*);

428: /* --------- options specific to particular eigensolvers -------- */

430: SLEPC_EXTERN PetscErrorCode NEPRIISetMaximumIterations(NEP,PetscInt);
431: SLEPC_EXTERN PetscErrorCode NEPRIIGetMaximumIterations(NEP,PetscInt*);
432: SLEPC_EXTERN PetscErrorCode NEPRIISetLagPreconditioner(NEP,PetscInt);
433: SLEPC_EXTERN PetscErrorCode NEPRIIGetLagPreconditioner(NEP,PetscInt*);
434: SLEPC_EXTERN PetscErrorCode NEPRIISetConstCorrectionTol(NEP,PetscBool);
435: SLEPC_EXTERN PetscErrorCode NEPRIIGetConstCorrectionTol(NEP,PetscBool*);
436: SLEPC_EXTERN PetscErrorCode NEPRIISetHermitian(NEP,PetscBool);
437: SLEPC_EXTERN PetscErrorCode NEPRIIGetHermitian(NEP,PetscBool*);
438: SLEPC_EXTERN PetscErrorCode NEPRIISetDeflationThreshold(NEP,PetscReal);
439: SLEPC_EXTERN PetscErrorCode NEPRIIGetDeflationThreshold(NEP,PetscReal*);
440: SLEPC_EXTERN PetscErrorCode NEPRIISetKSP(NEP,KSP);
441: SLEPC_EXTERN PetscErrorCode NEPRIIGetKSP(NEP,KSP*);

443: SLEPC_EXTERN PetscErrorCode NEPSLPSetDeflationThreshold(NEP,PetscReal);
444: SLEPC_EXTERN PetscErrorCode NEPSLPGetDeflationThreshold(NEP,PetscReal*);
445: SLEPC_EXTERN PetscErrorCode NEPSLPSetEPS(NEP,EPS);
446: SLEPC_EXTERN PetscErrorCode NEPSLPGetEPS(NEP,EPS*);
447: SLEPC_EXTERN PetscErrorCode NEPSLPSetEPSLeft(NEP,EPS);
448: SLEPC_EXTERN PetscErrorCode NEPSLPGetEPSLeft(NEP,EPS*);
449: SLEPC_EXTERN PetscErrorCode NEPSLPSetKSP(NEP,KSP);
450: SLEPC_EXTERN PetscErrorCode NEPSLPGetKSP(NEP,KSP*);

452: SLEPC_EXTERN PetscErrorCode NEPNArnoldiSetKSP(NEP,KSP);
453: SLEPC_EXTERN PetscErrorCode NEPNArnoldiGetKSP(NEP,KSP*);
454: SLEPC_EXTERN PetscErrorCode NEPNArnoldiSetLagPreconditioner(NEP,PetscInt);
455: SLEPC_EXTERN PetscErrorCode NEPNArnoldiGetLagPreconditioner(NEP,PetscInt*);

457: /*E
458:     NEPCISSExtraction - determines the extraction technique in the CISS solver

460:     Level: advanced

462: .seealso: NEPCISSSetExtraction(), NEPCISSGetExtraction()
463: E*/
464: typedef enum { NEP_CISS_EXTRACTION_RITZ,
465:                NEP_CISS_EXTRACTION_HANKEL,
466:                NEP_CISS_EXTRACTION_CAA    } NEPCISSExtraction;
467: SLEPC_EXTERN const char *NEPCISSExtractions[];

469: #if defined(PETSC_USE_COMPLEX) || defined(PETSC_CLANG_STATIC_ANALYZER)
470: SLEPC_EXTERN PetscErrorCode NEPCISSSetExtraction(NEP,NEPCISSExtraction);
471: SLEPC_EXTERN PetscErrorCode NEPCISSGetExtraction(NEP,NEPCISSExtraction*);
472: SLEPC_EXTERN PetscErrorCode NEPCISSSetSizes(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
473: SLEPC_EXTERN PetscErrorCode NEPCISSGetSizes(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
474: SLEPC_EXTERN PetscErrorCode NEPCISSSetThreshold(NEP,PetscReal,PetscReal);
475: SLEPC_EXTERN PetscErrorCode NEPCISSGetThreshold(NEP,PetscReal*,PetscReal*);
476: SLEPC_EXTERN PetscErrorCode NEPCISSSetRefinement(NEP,PetscInt,PetscInt);
477: SLEPC_EXTERN PetscErrorCode NEPCISSGetRefinement(NEP,PetscInt*,PetscInt*);
478: SLEPC_EXTERN PetscErrorCode NEPCISSGetKSPs(NEP,PetscInt*,KSP*[]);
479: #else
480: #define SlepcNEPCISSUnavailable(nep) do { \
481:     PetscFunctionBegin; \
482:     SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"%s() not available with real scalars",PETSC_FUNCTION_NAME); \
483:     } while (0)
484: static inline PetscErrorCode NEPCISSSetExtraction(NEP nep,PETSC_UNUSED NEPCISSExtraction ex) {SlepcNEPCISSUnavailable(nep);}
485: static inline PetscErrorCode NEPCISSGetExtraction(NEP nep,PETSC_UNUSED NEPCISSExtraction *ex) {SlepcNEPCISSUnavailable(nep);}
486: static inline PetscErrorCode NEPCISSSetSizes(NEP nep,PETSC_UNUSED PetscInt ip,PETSC_UNUSED PetscInt bs,PETSC_UNUSED PetscInt ms,PETSC_UNUSED PetscInt npart,PETSC_UNUSED PetscInt bsmax,PETSC_UNUSED PetscBool realmats) {SlepcNEPCISSUnavailable(nep);}
487: static inline PetscErrorCode NEPCISSGetSizes(NEP nep,PETSC_UNUSED PetscInt *ip,PETSC_UNUSED PetscInt *bs,PETSC_UNUSED PetscInt *ms,PETSC_UNUSED PetscInt *npart,PETSC_UNUSED PetscInt *bsmak,PETSC_UNUSED PetscBool *realmats) {SlepcNEPCISSUnavailable(nep);}
488: static inline PetscErrorCode NEPCISSSetThreshold(NEP nep,PETSC_UNUSED PetscReal delta,PETSC_UNUSED PetscReal spur) {SlepcNEPCISSUnavailable(nep);}
489: static inline PetscErrorCode NEPCISSGetThreshold(NEP nep,PETSC_UNUSED PetscReal *delta,PETSC_UNUSED PetscReal *spur) {SlepcNEPCISSUnavailable(nep);}
490: static inline PetscErrorCode NEPCISSSetRefinement(NEP nep,PETSC_UNUSED PetscInt inner,PETSC_UNUSED PetscInt blsize) {SlepcNEPCISSUnavailable(nep);}
491: static inline PetscErrorCode NEPCISSGetRefinement(NEP nep,PETSC_UNUSED PetscInt *inner,PETSC_UNUSED PetscInt *blsize) {SlepcNEPCISSUnavailable(nep);}
492: static inline PetscErrorCode NEPCISSGetKSPs(NEP nep,PETSC_UNUSED PetscInt *nsolve,PETSC_UNUSED KSP *ksp[]) {SlepcNEPCISSUnavailable(nep);}
493: #undef SlepcNEPCISSUnavailable
494: #endif

496: SLEPC_EXTERN PetscErrorCode NEPInterpolSetPEP(NEP,PEP);
497: SLEPC_EXTERN PetscErrorCode NEPInterpolGetPEP(NEP,PEP*);
498: SLEPC_EXTERN PetscErrorCode NEPInterpolSetInterpolation(NEP,PetscReal,PetscInt);
499: SLEPC_EXTERN PetscErrorCode NEPInterpolGetInterpolation(NEP,PetscReal*,PetscInt*);

501: /*S
502:   NEPNLEIGSSingularitiesFn - A prototype of a function that would be passed to NEPNLEIGSSetSingularitiesFunction()

504:   Calling Sequence:
505: +   nep   - eigensolver context obtained from NEPCreate()
506: .   maxnp - on input number of requested points in the discretization (can be set)
507: .   xi    - computed values of the discretization
508: -   ctx   - optional context, as set by NEPNLEIGSSetSingularitiesFunction()

510:   Level: advanced

512: .seealso: NEPNLEIGSSetSingularitiesFunction()
513: S*/
514: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode NEPNLEIGSSingularitiesFn(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *ctx);

516: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetSingularitiesFunction(NEP,NEPNLEIGSSingularitiesFn*,void*);
517: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetSingularitiesFunction(NEP,NEPNLEIGSSingularitiesFn**,void**);
518: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetRestart(NEP,PetscReal);
519: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetRestart(NEP,PetscReal*);
520: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetLocking(NEP,PetscBool);
521: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetLocking(NEP,PetscBool*);
522: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetInterpolation(NEP,PetscReal,PetscInt);
523: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetInterpolation(NEP,PetscReal*,PetscInt*);
524: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetRKShifts(NEP,PetscInt,PetscScalar[]);
525: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetRKShifts(NEP,PetscInt*,PetscScalar*[]);
526: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetKSPs(NEP,PetscInt*,KSP*[]);
527: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetFullBasis(NEP,PetscBool);
528: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetFullBasis(NEP,PetscBool*);
529: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetEPS(NEP,EPS);
530: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetEPS(NEP,EPS*);