Actual source code: slepcpep.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 polynomial eigenvalue solvers
 12: */

 14: #pragma once

 16: #include <slepceps.h>

 18: /* SUBMANSEC = PEP */

 20: SLEPC_EXTERN PetscErrorCode PEPInitializePackage(void);
 21: SLEPC_EXTERN PetscErrorCode PEPFinalizePackage(void);

 23: /*S
 24:      PEP - Abstract SLEPc object that manages all the polynomial eigenvalue
 25:      problem solvers.

 27:    Level: beginner

 29: .seealso:  `PEPCreate()`
 30: S*/
 31: typedef struct _p_PEP* PEP;

 33: /*J
 34:     PEPType - String with the name of a polynomial eigensolver

 36:    Level: beginner

 38: .seealso: `PEPSetType()`, `PEP`
 39: J*/
 40: typedef const char *PEPType;
 41: #define PEPLINEAR    "linear"
 42: #define PEPQARNOLDI  "qarnoldi"
 43: #define PEPTOAR      "toar"
 44: #define PEPSTOAR     "stoar"
 45: #define PEPJD        "jd"
 46: #define PEPCISS      "ciss"

 48: /* Logging support */
 49: SLEPC_EXTERN PetscClassId PEP_CLASSID;

 51: /*E
 52:     PEPProblemType - Determines the type of the polynomial eigenproblem

 54:     PEP_HERMITIAN is used when all A_i are Hermitian,
 55:     PEP_HYPERBOLIC is reserved for a QEP with Hermitian matrices, M>0, (x'Cx)^2 > 4(x'Mx)(x'Kx),
 56:     PEP_GYROSCOPIC is for aQEP with M, K  Hermitian, M>0, C skew-Hermitian.

 58:     Level: intermediate

 60: .seealso: `PEPSetProblemType()`, `PEPGetProblemType()`
 61: E*/
 62: typedef enum { PEP_GENERAL    = 1,
 63:                PEP_HERMITIAN  = 2,
 64:                PEP_HYPERBOLIC = 3,
 65:                PEP_GYROSCOPIC = 4
 66:              } PEPProblemType;

 68: /*E
 69:     PEPWhich - Determines which part of the spectrum is requested

 71:     Level: intermediate

 73: .seealso: `PEPSetWhichEigenpairs()`, `PEPGetWhichEigenpairs()`
 74: E*/
 75: typedef enum { PEP_LARGEST_MAGNITUDE  = 1,
 76:                PEP_SMALLEST_MAGNITUDE = 2,
 77:                PEP_LARGEST_REAL       = 3,
 78:                PEP_SMALLEST_REAL      = 4,
 79:                PEP_LARGEST_IMAGINARY  = 5,
 80:                PEP_SMALLEST_IMAGINARY = 6,
 81:                PEP_TARGET_MAGNITUDE   = 7,
 82:                PEP_TARGET_REAL        = 8,
 83:                PEP_TARGET_IMAGINARY   = 9,
 84:                PEP_ALL                = 10,
 85:                PEP_WHICH_USER         = 11 } PEPWhich;

 87: /*E
 88:     PEPBasis - The type of polynomial basis used to represent the polynomial
 89:     eigenproblem

 91:     Level: intermediate

 93: .seealso: `PEPSetBasis()`
 94: E*/
 95: typedef enum { PEP_BASIS_MONOMIAL,
 96:                PEP_BASIS_CHEBYSHEV1,
 97:                PEP_BASIS_CHEBYSHEV2,
 98:                PEP_BASIS_LEGENDRE,
 99:                PEP_BASIS_LAGUERRE,
100:                PEP_BASIS_HERMITE } PEPBasis;
101: SLEPC_EXTERN const char *PEPBasisTypes[];

103: /*E
104:     PEPScale - The scaling strategy

106:     Level: intermediate

108: .seealso: `PEPSetScale()`
109: E*/
110: typedef enum { PEP_SCALE_NONE,
111:                PEP_SCALE_SCALAR,
112:                PEP_SCALE_DIAGONAL,
113:                PEP_SCALE_BOTH } PEPScale;
114: SLEPC_EXTERN const char *PEPScaleTypes[];

116: /*E
117:     PEPRefine - The refinement type

119:     Level: intermediate

121: .seealso: `PEPSetRefine()`
122: E*/
123: typedef enum { PEP_REFINE_NONE,
124:                PEP_REFINE_SIMPLE,
125:                PEP_REFINE_MULTIPLE } PEPRefine;
126: SLEPC_EXTERN const char *PEPRefineTypes[];

128: /*E
129:     PEPRefineScheme - The scheme used for solving linear systems during iterative refinement

131:     Level: intermediate

133: .seealso: `PEPSetRefine()`
134: E*/
135: typedef enum { PEP_REFINE_SCHEME_SCHUR    = 1,
136:                PEP_REFINE_SCHEME_MBE      = 2,
137:                PEP_REFINE_SCHEME_EXPLICIT = 3 } PEPRefineScheme;
138: SLEPC_EXTERN const char *PEPRefineSchemes[];

140: /*E
141:     PEPExtract - The extraction type

143:     Level: intermediate

145: .seealso: `PEPSetExtract()`
146: E*/
147: typedef enum { PEP_EXTRACT_NONE       = 1,
148:                PEP_EXTRACT_NORM       = 2,
149:                PEP_EXTRACT_RESIDUAL   = 3,
150:                PEP_EXTRACT_STRUCTURED = 4 } PEPExtract;
151: SLEPC_EXTERN const char *PEPExtractTypes[];

153: /*E
154:     PEPErrorType - The error type used to assess accuracy of computed solutions

156:     Level: intermediate

158: .seealso: `PEPComputeError()`
159: E*/
160: typedef enum { PEP_ERROR_ABSOLUTE,
161:                PEP_ERROR_RELATIVE,
162:                PEP_ERROR_BACKWARD } PEPErrorType;
163: SLEPC_EXTERN const char *PEPErrorTypes[];

165: /*E
166:     PEPConv - Determines the convergence test

168:     Level: intermediate

170: .seealso: `PEPSetConvergenceTest()`, `PEPSetConvergenceTestFunction()`
171: E*/
172: typedef enum { PEP_CONV_ABS,
173:                PEP_CONV_REL,
174:                PEP_CONV_NORM,
175:                PEP_CONV_USER } PEPConv;

177: /*E
178:     PEPStop - Determines the stopping test

180:     Level: advanced

182: .seealso: `PEPSetStoppingTest()`, `PEPSetStoppingTestFunction()`
183: E*/
184: typedef enum { PEP_STOP_BASIC,
185:                PEP_STOP_USER } PEPStop;

187: /*E
188:     PEPConvergedReason - Reason an eigensolver was said to
189:          have converged or diverged

191:     Level: intermediate

193: .seealso: `PEPSolve()`, `PEPGetConvergedReason()`, `PEPSetTolerances()`
194: E*/
195: typedef enum {/* converged */
196:               PEP_CONVERGED_TOL                =  1,
197:               PEP_CONVERGED_USER               =  2,
198:               /* diverged */
199:               PEP_DIVERGED_ITS                 = -1,
200:               PEP_DIVERGED_BREAKDOWN           = -2,
201:               PEP_DIVERGED_SYMMETRY_LOST       = -3,
202:               PEP_CONVERGED_ITERATING          =  0} PEPConvergedReason;
203: SLEPC_EXTERN const char *const*PEPConvergedReasons;

205: SLEPC_EXTERN PetscErrorCode PEPCreate(MPI_Comm,PEP*);
206: SLEPC_EXTERN PetscErrorCode PEPDestroy(PEP*);
207: SLEPC_EXTERN PetscErrorCode PEPReset(PEP);
208: SLEPC_EXTERN PetscErrorCode PEPSetType(PEP,PEPType);
209: SLEPC_EXTERN PetscErrorCode PEPGetType(PEP,PEPType*);
210: SLEPC_EXTERN PetscErrorCode PEPSetProblemType(PEP,PEPProblemType);
211: SLEPC_EXTERN PetscErrorCode PEPGetProblemType(PEP,PEPProblemType*);
212: SLEPC_EXTERN PetscErrorCode PEPSetOperators(PEP,PetscInt,Mat[]);
213: SLEPC_EXTERN PetscErrorCode PEPGetOperators(PEP,PetscInt,Mat*);
214: SLEPC_EXTERN PetscErrorCode PEPGetNumMatrices(PEP,PetscInt*);
215: SLEPC_EXTERN PetscErrorCode PEPSetTarget(PEP,PetscScalar);
216: SLEPC_EXTERN PetscErrorCode PEPGetTarget(PEP,PetscScalar*);
217: SLEPC_EXTERN PetscErrorCode PEPSetInterval(PEP,PetscReal,PetscReal);
218: SLEPC_EXTERN PetscErrorCode PEPGetInterval(PEP,PetscReal*,PetscReal*);
219: SLEPC_EXTERN PetscErrorCode PEPSetFromOptions(PEP);
220: SLEPC_EXTERN PetscErrorCode PEPSetDSType(PEP);
221: SLEPC_EXTERN PetscErrorCode PEPSetUp(PEP);
222: SLEPC_EXTERN PetscErrorCode PEPSolve(PEP);
223: SLEPC_EXTERN PetscErrorCode PEPView(PEP,PetscViewer);
224: SLEPC_EXTERN PetscErrorCode PEPViewFromOptions(PEP,PetscObject,const char[]);
225: SLEPC_EXTERN PetscErrorCode PEPErrorView(PEP,PEPErrorType,PetscViewer);
226: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "PEPErrorView()", ) static inline PetscErrorCode PEPPrintSolution(PEP pep,PetscViewer v) {return PEPErrorView(pep,PEP_ERROR_BACKWARD,v);}
227: SLEPC_EXTERN PetscErrorCode PEPErrorViewFromOptions(PEP);
228: SLEPC_EXTERN PetscErrorCode PEPConvergedReasonView(PEP,PetscViewer);
229: SLEPC_EXTERN PetscErrorCode PEPConvergedReasonViewFromOptions(PEP);
230: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "PEPConvergedReasonView()", ) static inline PetscErrorCode PEPReasonView(PEP pep,PetscViewer v) {return PEPConvergedReasonView(pep,v);}
231: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "PEPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode PEPReasonViewFromOptions(PEP pep) {return PEPConvergedReasonViewFromOptions(pep);}
232: SLEPC_EXTERN PetscErrorCode PEPValuesView(PEP,PetscViewer);
233: SLEPC_EXTERN PetscErrorCode PEPValuesViewFromOptions(PEP);
234: SLEPC_EXTERN PetscErrorCode PEPVectorsView(PEP,PetscViewer);
235: SLEPC_EXTERN PetscErrorCode PEPVectorsViewFromOptions(PEP);
236: SLEPC_EXTERN PetscErrorCode PEPSetBV(PEP,BV);
237: SLEPC_EXTERN PetscErrorCode PEPGetBV(PEP,BV*);
238: SLEPC_EXTERN PetscErrorCode PEPSetRG(PEP,RG);
239: SLEPC_EXTERN PetscErrorCode PEPGetRG(PEP,RG*);
240: SLEPC_EXTERN PetscErrorCode PEPSetDS(PEP,DS);
241: SLEPC_EXTERN PetscErrorCode PEPGetDS(PEP,DS*);
242: SLEPC_EXTERN PetscErrorCode PEPSetST(PEP,ST);
243: SLEPC_EXTERN PetscErrorCode PEPGetST(PEP,ST*);
244: SLEPC_EXTERN PetscErrorCode PEPRefineGetKSP(PEP,KSP*);

246: SLEPC_EXTERN PetscErrorCode PEPSetTolerances(PEP,PetscReal,PetscInt);
247: SLEPC_EXTERN PetscErrorCode PEPGetTolerances(PEP,PetscReal*,PetscInt*);
248: SLEPC_EXTERN PetscErrorCode PEPGetConvergedReason(PEP,PEPConvergedReason*);

250: SLEPC_EXTERN PetscErrorCode PEPSetDimensions(PEP,PetscInt,PetscInt,PetscInt);
251: SLEPC_EXTERN PetscErrorCode PEPGetDimensions(PEP,PetscInt*,PetscInt*,PetscInt*);
252: SLEPC_EXTERN PetscErrorCode PEPSetScale(PEP,PEPScale,PetscReal,Vec,Vec,PetscInt,PetscReal);
253: SLEPC_EXTERN PetscErrorCode PEPGetScale(PEP,PEPScale*,PetscReal*,Vec*,Vec*,PetscInt*,PetscReal*);
254: SLEPC_EXTERN PetscErrorCode PEPSetRefine(PEP,PEPRefine,PetscInt,PetscReal,PetscInt,PEPRefineScheme);
255: SLEPC_EXTERN PetscErrorCode PEPGetRefine(PEP,PEPRefine*,PetscInt*,PetscReal*,PetscInt*,PEPRefineScheme*);
256: SLEPC_EXTERN PetscErrorCode PEPSetExtract(PEP,PEPExtract);
257: SLEPC_EXTERN PetscErrorCode PEPGetExtract(PEP,PEPExtract*);
258: SLEPC_EXTERN PetscErrorCode PEPSetBasis(PEP,PEPBasis);
259: SLEPC_EXTERN PetscErrorCode PEPGetBasis(PEP,PEPBasis*);

261: SLEPC_EXTERN PetscErrorCode PEPGetConverged(PEP,PetscInt*);
262: SLEPC_EXTERN PetscErrorCode PEPGetEigenpair(PEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
263: SLEPC_EXTERN PetscErrorCode PEPComputeError(PEP,PetscInt,PEPErrorType,PetscReal*);
264: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "PEPComputeError()", ) static inline PetscErrorCode PEPComputeRelativeError(PEP pep,PetscInt i,PetscReal *r) {return PEPComputeError(pep,i,PEP_ERROR_BACKWARD,r);}
265: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "PEPComputeError() with PEP_ERROR_ABSOLUTE", ) static inline PetscErrorCode PEPComputeResidualNorm(PEP pep,PetscInt i,PetscReal *r) {return PEPComputeError(pep,i,PEP_ERROR_ABSOLUTE,r);}
266: SLEPC_EXTERN PetscErrorCode PEPGetErrorEstimate(PEP,PetscInt,PetscReal*);
267: SLEPC_EXTERN PetscErrorCode PEPGetIterationNumber(PEP,PetscInt*);

269: SLEPC_EXTERN PetscErrorCode PEPSetInitialSpace(PEP,PetscInt,Vec[]);
270: SLEPC_EXTERN PetscErrorCode PEPSetWhichEigenpairs(PEP,PEPWhich);
271: SLEPC_EXTERN PetscErrorCode PEPGetWhichEigenpairs(PEP,PEPWhich*);

273: SLEPC_EXTERN PetscErrorCode PEPSetTrackAll(PEP,PetscBool);
274: SLEPC_EXTERN PetscErrorCode PEPGetTrackAll(PEP,PetscBool*);

276: /*S
277:   PEPMonitorFn - A function prototype for functions provided to PEPMonitorSet()

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

289:   Level: beginner

291: .seealso: `PEPMonitorSet()`
292: S*/
293: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PEPMonitorFn(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx);

295: /*S
296:   PEPMonitorRegisterFn - A function prototype for functions provided to PEPMonitorRegister()

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

308:   Level: beginner

310:   Note:
311:   This is an PEPMonitorFn specialized for a context of PetscViewerAndFormat.

313: .seealso: `PEPMonitorSet()`, `PEPMonitorRegister()`, `PEPMonitorFn`, `PEPMonitorRegisterCreateFn`, `PEPMonitorRegisterDestroyFn`
314: S*/
315: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PEPMonitorRegisterFn(PEP pep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,PetscViewerAndFormat *ctx);

317: /*S
318:   PEPMonitorRegisterCreateFn - A function prototype for functions that do the creation when provided to PEPMonitorRegister()

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

326:   Level: beginner

328: .seealso: `PEPMonitorRegisterFn`, `PEPMonitorSet()`, `PEPMonitorRegister()`, `PEPMonitorFn`, `PEPMonitorRegisterDestroyFn`
329: S*/
330: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PEPMonitorRegisterCreateFn(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **result);

332: /*S
333:   PEPMonitorRegisterDestroyFn - A function prototype for functions that do the after use destruction when provided to PEPMonitorRegister()

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

338:   Level: beginner

340: .seealso: `PEPMonitorRegisterFn`, `PEPMonitorSet()`, `PEPMonitorRegister()`, `PEPMonitorFn`, `PEPMonitorRegisterCreateFn`
341: S*/
342: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PEPMonitorRegisterDestroyFn(PetscViewerAndFormat **result);

344: SLEPC_EXTERN PetscErrorCode PEPMonitor(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
345: SLEPC_EXTERN PetscErrorCode PEPMonitorSet(PEP,PEPMonitorFn,void*,PetscCtxDestroyFn*);
346: SLEPC_EXTERN PetscErrorCode PEPMonitorCancel(PEP);
347: SLEPC_EXTERN PetscErrorCode PEPGetMonitorContext(PEP,void*);

349: SLEPC_EXTERN PetscErrorCode PEPMonitorSetFromOptions(PEP,const char[],const char[],void*,PetscBool);
350: SLEPC_EXTERN PEPMonitorRegisterFn        PEPMonitorFirst;
351: SLEPC_EXTERN PEPMonitorRegisterFn        PEPMonitorFirstDrawLG;
352: SLEPC_EXTERN PEPMonitorRegisterCreateFn  PEPMonitorFirstDrawLGCreate;
353: SLEPC_EXTERN PEPMonitorRegisterFn        PEPMonitorAll;
354: SLEPC_EXTERN PEPMonitorRegisterFn        PEPMonitorAllDrawLG;
355: SLEPC_EXTERN PEPMonitorRegisterCreateFn  PEPMonitorAllDrawLGCreate;
356: SLEPC_EXTERN PEPMonitorRegisterFn        PEPMonitorConverged;
357: SLEPC_EXTERN PEPMonitorRegisterCreateFn  PEPMonitorConvergedCreate;
358: SLEPC_EXTERN PEPMonitorRegisterFn        PEPMonitorConvergedDrawLG;
359: SLEPC_EXTERN PEPMonitorRegisterCreateFn  PEPMonitorConvergedDrawLGCreate;
360: SLEPC_EXTERN PEPMonitorRegisterDestroyFn PEPMonitorConvergedDestroy;

362: SLEPC_EXTERN PetscErrorCode PEPSetOptionsPrefix(PEP,const char*);
363: SLEPC_EXTERN PetscErrorCode PEPAppendOptionsPrefix(PEP,const char*);
364: SLEPC_EXTERN PetscErrorCode PEPGetOptionsPrefix(PEP,const char*[]);

366: SLEPC_EXTERN PetscFunctionList PEPList;
367: SLEPC_EXTERN PetscFunctionList PEPMonitorList;
368: SLEPC_EXTERN PetscFunctionList PEPMonitorCreateList;
369: SLEPC_EXTERN PetscFunctionList PEPMonitorDestroyList;
370: SLEPC_EXTERN PetscErrorCode PEPRegister(const char[],PetscErrorCode(*)(PEP));
371: SLEPC_EXTERN PetscErrorCode PEPMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PEPMonitorRegisterFn*,PEPMonitorRegisterCreateFn*,PEPMonitorRegisterDestroyFn*);

373: SLEPC_EXTERN PetscErrorCode PEPSetWorkVecs(PEP,PetscInt);
374: SLEPC_EXTERN PetscErrorCode PEPAllocateSolution(PEP,PetscInt);

376: /*S
377:   PEPConvergenceTestFn - A prototype of a PEP convergence test function that would be passed to PEPSetConvergenceTestFunction()

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

388:   Level: advanced

390: .seealso: `PEPSetConvergenceTestFunction()`
391: S*/
392: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PEPConvergenceTestFn(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx);

394: SLEPC_EXTERN PetscErrorCode PEPSetConvergenceTest(PEP,PEPConv);
395: SLEPC_EXTERN PetscErrorCode PEPGetConvergenceTest(PEP,PEPConv*);
396: SLEPC_EXTERN PEPConvergenceTestFn PEPConvergedAbsolute;
397: SLEPC_EXTERN PEPConvergenceTestFn PEPConvergedRelative;
398: SLEPC_EXTERN PEPConvergenceTestFn PEPConvergedNorm;
399: SLEPC_EXTERN PetscErrorCode PEPSetConvergenceTestFunction(PEP,PEPConvergenceTestFn*,void*,PetscCtxDestroyFn*);

401: /*S
402:   PEPStoppingTestFn - A prototype of a PEP stopping test function that would be passed to PEPSetStoppingTestFunction()

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

414:   Level: advanced

416: .seealso: `PEPSetStoppingTestFunction()`
417: S*/
418: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PEPStoppingTestFn(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx);

420: SLEPC_EXTERN PetscErrorCode PEPSetStoppingTest(PEP,PEPStop);
421: SLEPC_EXTERN PetscErrorCode PEPGetStoppingTest(PEP,PEPStop*);
422: SLEPC_EXTERN PEPStoppingTestFn PEPStoppingBasic;
423: SLEPC_EXTERN PetscErrorCode PEPSetStoppingTestFunction(PEP,PEPStoppingTestFn*,void*,PetscCtxDestroyFn*);

425: SLEPC_EXTERN PetscErrorCode PEPSetEigenvalueComparison(PEP,SlepcEigenvalueComparisonFn*,void*);

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

429: SLEPC_EXTERN PetscErrorCode PEPLinearSetLinearization(PEP,PetscReal,PetscReal);
430: SLEPC_EXTERN PetscErrorCode PEPLinearGetLinearization(PEP,PetscReal*,PetscReal*);
431: SLEPC_EXTERN PetscErrorCode PEPLinearSetExplicitMatrix(PEP,PetscBool);
432: SLEPC_EXTERN PetscErrorCode PEPLinearGetExplicitMatrix(PEP,PetscBool*);
433: SLEPC_EXTERN PetscErrorCode PEPLinearSetEPS(PEP,EPS);
434: SLEPC_EXTERN PetscErrorCode PEPLinearGetEPS(PEP,EPS*);
435: PETSC_DEPRECATED_FUNCTION(3, 10, 0, "PEPLinearSetLinearization()", ) static inline PetscErrorCode PEPLinearSetCompanionForm(PEP pep,PetscInt cform) {return (cform==1)?PEPLinearSetLinearization(pep,1.0,0.0):PEPLinearSetLinearization(pep,0.0,1.0);}
436: PETSC_DEPRECATED_FUNCTION(3, 10, 0, "PEPLinearGetLinearization()", ) static inline PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform) {(void)pep; if (cform) *cform=1; return PETSC_SUCCESS;}

438: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetRestart(PEP,PetscReal);
439: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetRestart(PEP,PetscReal*);
440: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetLocking(PEP,PetscBool);
441: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetLocking(PEP,PetscBool*);

443: SLEPC_EXTERN PetscErrorCode PEPTOARSetRestart(PEP,PetscReal);
444: SLEPC_EXTERN PetscErrorCode PEPTOARGetRestart(PEP,PetscReal*);
445: SLEPC_EXTERN PetscErrorCode PEPTOARSetLocking(PEP,PetscBool);
446: SLEPC_EXTERN PetscErrorCode PEPTOARGetLocking(PEP,PetscBool*);

448: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLinearization(PEP,PetscReal,PetscReal);
449: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLinearization(PEP,PetscReal*,PetscReal*);
450: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLocking(PEP,PetscBool);
451: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLocking(PEP,PetscBool*);
452: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDetectZeros(PEP,PetscBool);
453: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDetectZeros(PEP,PetscBool*);
454: SLEPC_EXTERN PetscErrorCode PEPSTOARGetInertias(PEP,PetscInt*,PetscReal*[],PetscInt*[]);
455: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDimensions(PEP,PetscInt,PetscInt,PetscInt);
456: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDimensions(PEP,PetscInt*,PetscInt*,PetscInt*);
457: SLEPC_EXTERN PetscErrorCode PEPSTOARSetCheckEigenvalueType(PEP,PetscBool);
458: SLEPC_EXTERN PetscErrorCode PEPSTOARGetCheckEigenvalueType(PEP,PetscBool*);
459: SLEPC_EXTERN PetscErrorCode PEPCheckDefiniteQEP(PEP,PetscReal*,PetscReal*,PetscInt*,PetscInt*);

461: /*E
462:     PEPJDProjection - The type of projection to be used in the Jacobi-Davidson solver

464:     Level: intermediate

466: .seealso: `PEPJDSetProjection()`
467: E*/
468: typedef enum { PEP_JD_PROJECTION_HARMONIC,
469:                PEP_JD_PROJECTION_ORTHOGONAL } PEPJDProjection;
470: SLEPC_EXTERN const char *PEPJDProjectionTypes[];

472: SLEPC_EXTERN PetscErrorCode PEPJDSetRestart(PEP,PetscReal);
473: SLEPC_EXTERN PetscErrorCode PEPJDGetRestart(PEP,PetscReal*);
474: SLEPC_EXTERN PetscErrorCode PEPJDSetFix(PEP,PetscReal);
475: SLEPC_EXTERN PetscErrorCode PEPJDGetFix(PEP,PetscReal*);
476: SLEPC_EXTERN PetscErrorCode PEPJDSetReusePreconditioner(PEP,PetscBool);
477: SLEPC_EXTERN PetscErrorCode PEPJDGetReusePreconditioner(PEP,PetscBool*);
478: SLEPC_EXTERN PetscErrorCode PEPJDSetMinimalityIndex(PEP,PetscInt);
479: SLEPC_EXTERN PetscErrorCode PEPJDGetMinimalityIndex(PEP,PetscInt*);
480: SLEPC_EXTERN PetscErrorCode PEPJDSetProjection(PEP,PEPJDProjection);
481: SLEPC_EXTERN PetscErrorCode PEPJDGetProjection(PEP,PEPJDProjection*);

483: /*E
484:     PEPCISSExtraction - determines the extraction technique in the CISS solver

486:     Level: advanced

488: .seealso: `PEPCISSSetExtraction()`, `PEPCISSGetExtraction()`
489: E*/
490: typedef enum { PEP_CISS_EXTRACTION_RITZ,
491:                PEP_CISS_EXTRACTION_HANKEL,
492:                PEP_CISS_EXTRACTION_CAA    } PEPCISSExtraction;
493: SLEPC_EXTERN const char *PEPCISSExtractions[];

495: #if defined(PETSC_USE_COMPLEX) || defined(PETSC_CLANG_STATIC_ANALYZER)
496: SLEPC_EXTERN PetscErrorCode PEPCISSSetExtraction(PEP,PEPCISSExtraction);
497: SLEPC_EXTERN PetscErrorCode PEPCISSGetExtraction(PEP,PEPCISSExtraction*);
498: SLEPC_EXTERN PetscErrorCode PEPCISSSetSizes(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
499: SLEPC_EXTERN PetscErrorCode PEPCISSGetSizes(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
500: SLEPC_EXTERN PetscErrorCode PEPCISSSetThreshold(PEP,PetscReal,PetscReal);
501: SLEPC_EXTERN PetscErrorCode PEPCISSGetThreshold(PEP,PetscReal*,PetscReal*);
502: SLEPC_EXTERN PetscErrorCode PEPCISSSetRefinement(PEP,PetscInt,PetscInt);
503: SLEPC_EXTERN PetscErrorCode PEPCISSGetRefinement(PEP,PetscInt*,PetscInt*);
504: SLEPC_EXTERN PetscErrorCode PEPCISSGetKSPs(PEP,PetscInt*,KSP*[]);
505: #else
506: #define SlepcPEPCISSUnavailable(pep) do { \
507:     PetscFunctionBegin; \
508:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"%s() not available with real scalars",PETSC_FUNCTION_NAME); \
509:     } while (0)
510: static inline PetscErrorCode PEPCISSSetExtraction(PEP pep,PETSC_UNUSED PEPCISSExtraction ex) {SlepcPEPCISSUnavailable(pep);}
511: static inline PetscErrorCode PEPCISSGetExtraction(PEP pep,PETSC_UNUSED PEPCISSExtraction *ex) {SlepcPEPCISSUnavailable(pep);}
512: static inline PetscErrorCode PEPCISSSetSizes(PEP pep,PETSC_UNUSED PetscInt ip,PETSC_UNUSED PetscInt bs,PETSC_UNUSED PetscInt ms,PETSC_UNUSED PetscInt npart,PETSC_UNUSED PetscInt bsmax,PETSC_UNUSED PetscBool realmats) {SlepcPEPCISSUnavailable(pep);}
513: static inline PetscErrorCode PEPCISSGetSizes(PEP pep,PETSC_UNUSED PetscInt *ip,PETSC_UNUSED PetscInt *bs,PETSC_UNUSED PetscInt *ms,PETSC_UNUSED PetscInt *npart,PETSC_UNUSED PetscInt *bsmak,PETSC_UNUSED PetscBool *realmats) {SlepcPEPCISSUnavailable(pep);}
514: static inline PetscErrorCode PEPCISSSetThreshold(PEP pep,PETSC_UNUSED PetscReal delta,PETSC_UNUSED PetscReal spur) {SlepcPEPCISSUnavailable(pep);}
515: static inline PetscErrorCode PEPCISSGetThreshold(PEP pep,PETSC_UNUSED PetscReal *delta,PETSC_UNUSED PetscReal *spur) {SlepcPEPCISSUnavailable(pep);}
516: static inline PetscErrorCode PEPCISSSetRefinement(PEP pep,PETSC_UNUSED PetscInt inner,PETSC_UNUSED PetscInt blsize) {SlepcPEPCISSUnavailable(pep);}
517: static inline PetscErrorCode PEPCISSGetRefinement(PEP pep,PETSC_UNUSED PetscInt *inner,PETSC_UNUSED PetscInt *blsize) {SlepcPEPCISSUnavailable(pep);}
518: static inline PetscErrorCode PEPCISSGetKSPs(PEP pep,PETSC_UNUSED PetscInt *nsolve,PETSC_UNUSED KSP *ksp[]) {SlepcPEPCISSUnavailable(pep);}
519: #undef SlepcPEPCISSUnavailable
520: #endif