Actual source code: slepcpep.h

slepc-main 2024-11-09
Report Typos and Errors
  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:     Level: intermediate

 56: .seealso: PEPSetProblemType(), PEPGetProblemType()
 57: E*/
 58: typedef enum { PEP_GENERAL=1,
 59:                PEP_HERMITIAN,   /* All A_i  Hermitian */
 60:                PEP_HYPERBOLIC,  /* QEP with Hermitian matrices, M>0, (x'Cx)^2 > 4(x'Mx)(x'Kx) */
 61:                PEP_GYROSCOPIC   /* QEP with M, K  Hermitian, M>0, C skew-Hermitian */
 62:              } PEPProblemType;

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

 67:     Level: intermediate

 69: .seealso: PEPSetWhichEigenpairs(), PEPGetWhichEigenpairs()
 70: E*/
 71: typedef enum { PEP_LARGEST_MAGNITUDE=1,
 72:                PEP_SMALLEST_MAGNITUDE,
 73:                PEP_LARGEST_REAL,
 74:                PEP_SMALLEST_REAL,
 75:                PEP_LARGEST_IMAGINARY,
 76:                PEP_SMALLEST_IMAGINARY,
 77:                PEP_TARGET_MAGNITUDE,
 78:                PEP_TARGET_REAL,
 79:                PEP_TARGET_IMAGINARY,
 80:                PEP_ALL,
 81:                PEP_WHICH_USER } PEPWhich;

 83: /*E
 84:     PEPBasis - The type of polynomial basis used to represent the polynomial
 85:     eigenproblem

 87:     Level: intermediate

 89: .seealso: PEPSetBasis()
 90: E*/
 91: typedef enum { PEP_BASIS_MONOMIAL,
 92:                PEP_BASIS_CHEBYSHEV1,
 93:                PEP_BASIS_CHEBYSHEV2,
 94:                PEP_BASIS_LEGENDRE,
 95:                PEP_BASIS_LAGUERRE,
 96:                PEP_BASIS_HERMITE } PEPBasis;
 97: SLEPC_EXTERN const char *PEPBasisTypes[];

 99: /*E
100:     PEPScale - The scaling strategy

102:     Level: intermediate

104: .seealso: PEPSetScale()
105: E*/
106: typedef enum { PEP_SCALE_NONE,
107:                PEP_SCALE_SCALAR,
108:                PEP_SCALE_DIAGONAL,
109:                PEP_SCALE_BOTH } PEPScale;
110: SLEPC_EXTERN const char *PEPScaleTypes[];

112: /*E
113:     PEPRefine - The refinement type

115:     Level: intermediate

117: .seealso: PEPSetRefine()
118: E*/
119: typedef enum { PEP_REFINE_NONE,
120:                PEP_REFINE_SIMPLE,
121:                PEP_REFINE_MULTIPLE } PEPRefine;
122: SLEPC_EXTERN const char *PEPRefineTypes[];

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

127:     Level: intermediate

129: .seealso: PEPSetRefine()
130: E*/
131: typedef enum { PEP_REFINE_SCHEME_SCHUR=1,
132:                PEP_REFINE_SCHEME_MBE,
133:                PEP_REFINE_SCHEME_EXPLICIT } PEPRefineScheme;
134: SLEPC_EXTERN const char *PEPRefineSchemes[];

136: /*E
137:     PEPExtract - The extraction type

139:     Level: intermediate

141: .seealso: PEPSetExtract()
142: E*/
143: typedef enum { PEP_EXTRACT_NONE=1,
144:                PEP_EXTRACT_NORM,
145:                PEP_EXTRACT_RESIDUAL,
146:                PEP_EXTRACT_STRUCTURED } PEPExtract;
147: SLEPC_EXTERN const char *PEPExtractTypes[];

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

152:     Level: intermediate

154: .seealso: PEPComputeError()
155: E*/
156: typedef enum { PEP_ERROR_ABSOLUTE,
157:                PEP_ERROR_RELATIVE,
158:                PEP_ERROR_BACKWARD } PEPErrorType;
159: SLEPC_EXTERN const char *PEPErrorTypes[];

161: /*E
162:     PEPConv - Determines the convergence test

164:     Level: intermediate

166: .seealso: PEPSetConvergenceTest(), PEPSetConvergenceTestFunction()
167: E*/
168: typedef enum { PEP_CONV_ABS,
169:                PEP_CONV_REL,
170:                PEP_CONV_NORM,
171:                PEP_CONV_USER } PEPConv;

173: /*E
174:     PEPStop - Determines the stopping test

176:     Level: advanced

178: .seealso: PEPSetStoppingTest(), PEPSetStoppingTestFunction()
179: E*/
180: typedef enum { PEP_STOP_BASIC,
181:                PEP_STOP_USER } PEPStop;

183: /*E
184:     PEPConvergedReason - Reason an eigensolver was said to
185:          have converged or diverged

187:     Level: intermediate

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

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

242: SLEPC_EXTERN PetscErrorCode PEPSetTolerances(PEP,PetscReal,PetscInt);
243: SLEPC_EXTERN PetscErrorCode PEPGetTolerances(PEP,PetscReal*,PetscInt*);
244: SLEPC_EXTERN PetscErrorCode PEPSetConvergenceTest(PEP,PEPConv);
245: SLEPC_EXTERN PetscErrorCode PEPGetConvergenceTest(PEP,PEPConv*);
246: SLEPC_EXTERN PetscErrorCode PEPConvergedAbsolute(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
247: SLEPC_EXTERN PetscErrorCode PEPConvergedRelative(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
248: SLEPC_EXTERN PetscErrorCode PEPConvergedNorm(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
249: SLEPC_EXTERN PetscErrorCode PEPSetStoppingTest(PEP,PEPStop);
250: SLEPC_EXTERN PetscErrorCode PEPGetStoppingTest(PEP,PEPStop*);
251: SLEPC_EXTERN PetscErrorCode PEPStoppingBasic(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*);
252: SLEPC_EXTERN PetscErrorCode PEPGetConvergedReason(PEP,PEPConvergedReason*);

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

265: SLEPC_EXTERN PetscErrorCode PEPGetConverged(PEP,PetscInt*);
266: SLEPC_EXTERN PetscErrorCode PEPGetEigenpair(PEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
267: SLEPC_EXTERN PetscErrorCode PEPComputeError(PEP,PetscInt,PEPErrorType,PetscReal*);
268: 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);}
269: 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);}
270: SLEPC_EXTERN PetscErrorCode PEPGetErrorEstimate(PEP,PetscInt,PetscReal*);
271: SLEPC_EXTERN PetscErrorCode PEPGetIterationNumber(PEP,PetscInt*);

273: SLEPC_EXTERN PetscErrorCode PEPSetInitialSpace(PEP,PetscInt,Vec[]);
274: SLEPC_EXTERN PetscErrorCode PEPSetWhichEigenpairs(PEP,PEPWhich);
275: SLEPC_EXTERN PetscErrorCode PEPGetWhichEigenpairs(PEP,PEPWhich*);

277: SLEPC_EXTERN PetscErrorCode PEPSetTrackAll(PEP,PetscBool);
278: SLEPC_EXTERN PetscErrorCode PEPGetTrackAll(PEP,PetscBool*);

280: SLEPC_EXTERN PetscErrorCode PEPMonitor(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
281: SLEPC_EXTERN PetscErrorCode PEPMonitorSet(PEP,PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscCtxDestroyFn*);
282: SLEPC_EXTERN PetscErrorCode PEPMonitorCancel(PEP);
283: SLEPC_EXTERN PetscErrorCode PEPGetMonitorContext(PEP,void*);

285: SLEPC_EXTERN PetscErrorCode PEPMonitorSetFromOptions(PEP,const char[],const char[],void*,PetscBool);
286: SLEPC_EXTERN PetscErrorCode PEPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char*[],int,int,int,int,PetscDrawLG*);
287: SLEPC_EXTERN PetscErrorCode PEPMonitorFirst(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
288: SLEPC_EXTERN PetscErrorCode PEPMonitorFirstDrawLG(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
289: SLEPC_EXTERN PetscErrorCode PEPMonitorFirstDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
290: SLEPC_EXTERN PetscErrorCode PEPMonitorAll(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
291: SLEPC_EXTERN PetscErrorCode PEPMonitorAllDrawLG(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
292: SLEPC_EXTERN PetscErrorCode PEPMonitorAllDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
293: SLEPC_EXTERN PetscErrorCode PEPMonitorConverged(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
294: SLEPC_EXTERN PetscErrorCode PEPMonitorConvergedCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
295: SLEPC_EXTERN PetscErrorCode PEPMonitorConvergedDrawLG(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
296: SLEPC_EXTERN PetscErrorCode PEPMonitorConvergedDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
297: SLEPC_EXTERN PetscErrorCode PEPMonitorConvergedDestroy(PetscViewerAndFormat**);

299: SLEPC_EXTERN PetscErrorCode PEPSetOptionsPrefix(PEP,const char*);
300: SLEPC_EXTERN PetscErrorCode PEPAppendOptionsPrefix(PEP,const char*);
301: SLEPC_EXTERN PetscErrorCode PEPGetOptionsPrefix(PEP,const char*[]);

303: SLEPC_EXTERN PetscFunctionList PEPList;
304: SLEPC_EXTERN PetscFunctionList PEPMonitorList;
305: SLEPC_EXTERN PetscFunctionList PEPMonitorCreateList;
306: SLEPC_EXTERN PetscFunctionList PEPMonitorDestroyList;
307: SLEPC_EXTERN PetscErrorCode PEPRegister(const char[],PetscErrorCode(*)(PEP));
308: SLEPC_EXTERN PetscErrorCode PEPMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));

310: SLEPC_EXTERN PetscErrorCode PEPSetWorkVecs(PEP,PetscInt);
311: SLEPC_EXTERN PetscErrorCode PEPAllocateSolution(PEP,PetscInt);

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

316:   Calling Sequence:
317: +   pep    - eigensolver context obtained from PEPCreate()
318: .   eigr   - real part of the eigenvalue
319: .   eigi   - imaginary part of the eigenvalue
320: .   res    - residual norm associated to the eigenpair
321: .   errest - [output] computed error estimate
322: -   ctx    - [optional] user-defined context for private data for the
323:              convergence test routine (may be NULL)

325:   Level: advanced

327: .seealso: PEPSetConvergenceTestFunction()
328: S*/
329: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(PEPConvergenceTestFn)(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx);

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

334:   Calling Sequence:
335: +   pep    - eigensolver context obtained from PEPCreate()
336: .   its    - current number of iterations
337: .   max_it - maximum number of iterations
338: .   nconv  - number of currently converged eigenpairs
339: .   nev    - number of requested eigenpairs
340: .   reason - [output] result of the stopping test
341: -   ctx    - [optional] user-defined context for private data for the
342:              stopping test routine (may be NULL)

344:   Level: advanced

346: .seealso: PEPSetStoppingTestFunction()
347: S*/
348: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(PEPStoppingTestFn)(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx);

350: SLEPC_EXTERN PetscErrorCode PEPSetConvergenceTestFunction(PEP,PEPConvergenceTestFn*,void*,PetscCtxDestroyFn*);
351: SLEPC_EXTERN PetscErrorCode PEPSetStoppingTestFunction(PEP,PEPStoppingTestFn*,void*,PetscCtxDestroyFn*);
352: SLEPC_EXTERN PetscErrorCode PEPSetEigenvalueComparison(PEP,SlepcEigenvalueComparisonFn*,void*);

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

356: SLEPC_EXTERN PetscErrorCode PEPLinearSetLinearization(PEP,PetscReal,PetscReal);
357: SLEPC_EXTERN PetscErrorCode PEPLinearGetLinearization(PEP,PetscReal*,PetscReal*);
358: SLEPC_EXTERN PetscErrorCode PEPLinearSetExplicitMatrix(PEP,PetscBool);
359: SLEPC_EXTERN PetscErrorCode PEPLinearGetExplicitMatrix(PEP,PetscBool*);
360: SLEPC_EXTERN PetscErrorCode PEPLinearSetEPS(PEP,EPS);
361: SLEPC_EXTERN PetscErrorCode PEPLinearGetEPS(PEP,EPS*);
362: 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);}
363: PETSC_DEPRECATED_FUNCTION(3, 10, 0, "PEPLinearGetLinearization()", ) static inline PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform) {(void)pep; if (cform) *cform=1; return PETSC_SUCCESS;}

365: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetRestart(PEP,PetscReal);
366: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetRestart(PEP,PetscReal*);
367: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetLocking(PEP,PetscBool);
368: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetLocking(PEP,PetscBool*);

370: SLEPC_EXTERN PetscErrorCode PEPTOARSetRestart(PEP,PetscReal);
371: SLEPC_EXTERN PetscErrorCode PEPTOARGetRestart(PEP,PetscReal*);
372: SLEPC_EXTERN PetscErrorCode PEPTOARSetLocking(PEP,PetscBool);
373: SLEPC_EXTERN PetscErrorCode PEPTOARGetLocking(PEP,PetscBool*);

375: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLinearization(PEP,PetscReal,PetscReal);
376: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLinearization(PEP,PetscReal*,PetscReal*);
377: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLocking(PEP,PetscBool);
378: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLocking(PEP,PetscBool*);
379: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDetectZeros(PEP,PetscBool);
380: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDetectZeros(PEP,PetscBool*);
381: SLEPC_EXTERN PetscErrorCode PEPSTOARGetInertias(PEP,PetscInt*,PetscReal**,PetscInt**);
382: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDimensions(PEP,PetscInt,PetscInt,PetscInt);
383: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDimensions(PEP,PetscInt*,PetscInt*,PetscInt*);
384: SLEPC_EXTERN PetscErrorCode PEPSTOARSetCheckEigenvalueType(PEP,PetscBool);
385: SLEPC_EXTERN PetscErrorCode PEPSTOARGetCheckEigenvalueType(PEP,PetscBool*);
386: SLEPC_EXTERN PetscErrorCode PEPCheckDefiniteQEP(PEP,PetscReal*,PetscReal*,PetscInt*,PetscInt*);

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

391:     Level: intermediate

393: .seealso: PEPJDSetProjection()
394: E*/
395: typedef enum { PEP_JD_PROJECTION_HARMONIC,
396:                PEP_JD_PROJECTION_ORTHOGONAL } PEPJDProjection;
397: SLEPC_EXTERN const char *PEPJDProjectionTypes[];

399: SLEPC_EXTERN PetscErrorCode PEPJDSetRestart(PEP,PetscReal);
400: SLEPC_EXTERN PetscErrorCode PEPJDGetRestart(PEP,PetscReal*);
401: SLEPC_EXTERN PetscErrorCode PEPJDSetFix(PEP,PetscReal);
402: SLEPC_EXTERN PetscErrorCode PEPJDGetFix(PEP,PetscReal*);
403: SLEPC_EXTERN PetscErrorCode PEPJDSetReusePreconditioner(PEP,PetscBool);
404: SLEPC_EXTERN PetscErrorCode PEPJDGetReusePreconditioner(PEP,PetscBool*);
405: SLEPC_EXTERN PetscErrorCode PEPJDSetMinimalityIndex(PEP,PetscInt);
406: SLEPC_EXTERN PetscErrorCode PEPJDGetMinimalityIndex(PEP,PetscInt*);
407: SLEPC_EXTERN PetscErrorCode PEPJDSetProjection(PEP,PEPJDProjection);
408: SLEPC_EXTERN PetscErrorCode PEPJDGetProjection(PEP,PEPJDProjection*);

410: /*E
411:     PEPCISSExtraction - determines the extraction technique in the CISS solver

413:     Level: advanced

415: .seealso: PEPCISSSetExtraction(), PEPCISSGetExtraction()
416: E*/
417: typedef enum { PEP_CISS_EXTRACTION_RITZ,
418:                PEP_CISS_EXTRACTION_HANKEL,
419:                PEP_CISS_EXTRACTION_CAA    } PEPCISSExtraction;
420: SLEPC_EXTERN const char *PEPCISSExtractions[];

422: #if defined(PETSC_USE_COMPLEX) || defined(PETSC_CLANG_STATIC_ANALYZER)
423: SLEPC_EXTERN PetscErrorCode PEPCISSSetExtraction(PEP,PEPCISSExtraction);
424: SLEPC_EXTERN PetscErrorCode PEPCISSGetExtraction(PEP,PEPCISSExtraction*);
425: SLEPC_EXTERN PetscErrorCode PEPCISSSetSizes(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
426: SLEPC_EXTERN PetscErrorCode PEPCISSGetSizes(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
427: SLEPC_EXTERN PetscErrorCode PEPCISSSetThreshold(PEP,PetscReal,PetscReal);
428: SLEPC_EXTERN PetscErrorCode PEPCISSGetThreshold(PEP,PetscReal*,PetscReal*);
429: SLEPC_EXTERN PetscErrorCode PEPCISSSetRefinement(PEP,PetscInt,PetscInt);
430: SLEPC_EXTERN PetscErrorCode PEPCISSGetRefinement(PEP,PetscInt*,PetscInt*);
431: SLEPC_EXTERN PetscErrorCode PEPCISSGetKSPs(PEP,PetscInt*,KSP**);
432: #else
433: #define SlepcPEPCISSUnavailable(pep) do { \
434:     PetscFunctionBegin; \
435:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"%s() not available with real scalars",PETSC_FUNCTION_NAME); \
436:     } while (0)
437: static inline PetscErrorCode PEPCISSSetExtraction(PEP pep,PETSC_UNUSED PEPCISSExtraction ex) {SlepcPEPCISSUnavailable(pep);}
438: static inline PetscErrorCode PEPCISSGetExtraction(PEP pep,PETSC_UNUSED PEPCISSExtraction *ex) {SlepcPEPCISSUnavailable(pep);}
439: 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);}
440: 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);}
441: static inline PetscErrorCode PEPCISSSetThreshold(PEP pep,PETSC_UNUSED PetscReal delta,PETSC_UNUSED PetscReal spur) {SlepcPEPCISSUnavailable(pep);}
442: static inline PetscErrorCode PEPCISSGetThreshold(PEP pep,PETSC_UNUSED PetscReal *delta,PETSC_UNUSED PetscReal *spur) {SlepcPEPCISSUnavailable(pep);}
443: static inline PetscErrorCode PEPCISSSetRefinement(PEP pep,PETSC_UNUSED PetscInt inner,PETSC_UNUSED PetscInt blsize) {SlepcPEPCISSUnavailable(pep);}
444: static inline PetscErrorCode PEPCISSGetRefinement(PEP pep,PETSC_UNUSED PetscInt *inner,PETSC_UNUSED PetscInt *blsize) {SlepcPEPCISSUnavailable(pep);}
445: static inline PetscErrorCode PEPCISSGetKSPs(PEP pep,PETSC_UNUSED PetscInt *nsolve,PETSC_UNUSED KSP **ksp) {SlepcPEPCISSUnavailable(pep);}
446: #undef SlepcPEPCISSUnavailable
447: #endif