Actual source code: slepcpep.h

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: #if !defined(SLEPCPEP_H)
 15: #define SLEPCPEP_H
 16: #include <slepceps.h>

 18: SLEPC_EXTERN PetscErrorCode PEPInitializePackage(void);

 20: /*S
 21:      PEP - Abstract SLEPc object that manages all the polynomial eigenvalue
 22:      problem solvers.

 24:    Level: beginner

 26: .seealso:  PEPCreate()
 27: S*/
 28: typedef struct _p_PEP* PEP;

 30: /*J
 31:     PEPType - String with the name of a polynomial eigensolver

 33:    Level: beginner

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

 45: /* Logging support */
 46: SLEPC_EXTERN PetscClassId PEP_CLASSID;

 48: /*E
 49:     PEPProblemType - Determines the type of the polynomial eigenproblem

 51:     Level: intermediate

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

 61: /*E
 62:     PEPWhich - Determines which part of the spectrum is requested

 64:     Level: intermediate

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

 80: /*E
 81:     PEPBasis - The type of polynomial basis used to represent the polynomial
 82:     eigenproblem

 84:     Level: intermediate

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

 96: /*E
 97:     PEPScale - The scaling strategy

 99:     Level: intermediate

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

109: /*E
110:     PEPRefine - The refinement type

112:     Level: intermediate

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

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

124:     Level: intermediate

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

133: /*E
134:     PEPExtract - The extraction type

136:     Level: intermediate

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

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

149:     Level: intermediate

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

158: /*E
159:     PEPConv - Determines the convergence test

161:     Level: intermediate

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

170: /*E
171:     PEPStop - Determines the stopping test

173:     Level: advanced

175: .seealso: PEPSetStoppingTest(), PEPSetStoppingTestFunction()
176: E*/
177: typedef enum { PEP_STOP_BASIC,
178:                PEP_STOP_USER } PEPStop;

180: /*E
181:     PEPConvergedReason - Reason an eigensolver was said to
182:          have converged or diverged

184:     Level: intermediate

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

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

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

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

263: SLEPC_EXTERN PetscErrorCode PEPGetConverged(PEP,PetscInt*);
264: SLEPC_EXTERN PetscErrorCode PEPGetEigenpair(PEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
265: SLEPC_EXTERN PetscErrorCode PEPComputeError(PEP,PetscInt,PEPErrorType,PetscReal*);
266: PETSC_DEPRECATED_FUNCTION("Use PEPComputeError()") PETSC_STATIC_INLINE PetscErrorCode PEPComputeRelativeError(PEP pep,PetscInt i,PetscReal *r) {return PEPComputeError(pep,i,PEP_ERROR_BACKWARD,r);}
267: PETSC_DEPRECATED_FUNCTION("Use PEPComputeError() with PEP_ERROR_ABSOLUTE") PETSC_STATIC_INLINE PetscErrorCode PEPComputeResidualNorm(PEP pep,PetscInt i,PetscReal *r) {return PEPComputeError(pep,i,PEP_ERROR_ABSOLUTE,r);}
268: SLEPC_EXTERN PetscErrorCode PEPGetErrorEstimate(PEP,PetscInt,PetscReal*);
269: SLEPC_EXTERN PetscErrorCode PEPGetIterationNumber(PEP,PetscInt*);

271: SLEPC_EXTERN PetscErrorCode PEPSetInitialSpace(PEP,PetscInt,Vec[]);
272: SLEPC_EXTERN PetscErrorCode PEPSetWhichEigenpairs(PEP,PEPWhich);
273: SLEPC_EXTERN PetscErrorCode PEPGetWhichEigenpairs(PEP,PEPWhich*);
274: SLEPC_EXTERN PetscErrorCode PEPSetEigenvalueComparison(PEP,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);

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

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

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

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

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

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

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

314: SLEPC_EXTERN PetscErrorCode PEPLinearSetLinearization(PEP,PetscReal,PetscReal);
315: SLEPC_EXTERN PetscErrorCode PEPLinearGetLinearization(PEP,PetscReal*,PetscReal*);
316: SLEPC_EXTERN PetscErrorCode PEPLinearSetExplicitMatrix(PEP,PetscBool);
317: SLEPC_EXTERN PetscErrorCode PEPLinearGetExplicitMatrix(PEP,PetscBool*);
318: SLEPC_EXTERN PetscErrorCode PEPLinearSetEPS(PEP,EPS);
319: SLEPC_EXTERN PetscErrorCode PEPLinearGetEPS(PEP,EPS*);
320: PETSC_DEPRECATED_FUNCTION("Use PEPLinearSetLinearization()") PETSC_STATIC_INLINE PetscErrorCode PEPLinearSetCompanionForm(PEP pep,PetscInt cform) {return (cform==1)?PEPLinearSetLinearization(pep,1.0,0.0):PEPLinearSetLinearization(pep,0.0,1.0);}
321: PETSC_DEPRECATED_FUNCTION("Use PEPLinearGetLinearization()") PETSC_STATIC_INLINE PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform) {(void)pep; if (cform) *cform=1; return 0;}

323: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetRestart(PEP,PetscReal);
324: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetRestart(PEP,PetscReal*);
325: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetLocking(PEP,PetscBool);
326: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetLocking(PEP,PetscBool*);

328: SLEPC_EXTERN PetscErrorCode PEPTOARSetRestart(PEP,PetscReal);
329: SLEPC_EXTERN PetscErrorCode PEPTOARGetRestart(PEP,PetscReal*);
330: SLEPC_EXTERN PetscErrorCode PEPTOARSetLocking(PEP,PetscBool);
331: SLEPC_EXTERN PetscErrorCode PEPTOARGetLocking(PEP,PetscBool*);

333: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLinearization(PEP,PetscReal,PetscReal);
334: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLinearization(PEP,PetscReal*,PetscReal*);
335: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLocking(PEP,PetscBool);
336: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLocking(PEP,PetscBool*);
337: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDetectZeros(PEP,PetscBool);
338: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDetectZeros(PEP,PetscBool*);
339: SLEPC_EXTERN PetscErrorCode PEPSTOARGetInertias(PEP,PetscInt*,PetscReal**,PetscInt**);
340: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDimensions(PEP,PetscInt,PetscInt,PetscInt);
341: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDimensions(PEP,PetscInt*,PetscInt*,PetscInt*);
342: SLEPC_EXTERN PetscErrorCode PEPSTOARSetCheckEigenvalueType(PEP,PetscBool);
343: SLEPC_EXTERN PetscErrorCode PEPSTOARGetCheckEigenvalueType(PEP,PetscBool*);
344: SLEPC_EXTERN PetscErrorCode PEPCheckDefiniteQEP(PEP,PetscReal*,PetscReal*,PetscInt*,PetscInt*);

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

349:     Level: intermediate

351: .seealso: PEPJDSetProjection()
352: E*/
353: typedef enum { PEP_JD_PROJECTION_HARMONIC,
354:                PEP_JD_PROJECTION_ORTHOGONAL } PEPJDProjection;
355: SLEPC_EXTERN const char *PEPJDProjectionTypes[];

357: SLEPC_EXTERN PetscErrorCode PEPJDSetRestart(PEP,PetscReal);
358: SLEPC_EXTERN PetscErrorCode PEPJDGetRestart(PEP,PetscReal*);
359: SLEPC_EXTERN PetscErrorCode PEPJDSetFix(PEP,PetscReal);
360: SLEPC_EXTERN PetscErrorCode PEPJDGetFix(PEP,PetscReal*);
361: SLEPC_EXTERN PetscErrorCode PEPJDSetReusePreconditioner(PEP,PetscBool);
362: SLEPC_EXTERN PetscErrorCode PEPJDGetReusePreconditioner(PEP,PetscBool*);
363: SLEPC_EXTERN PetscErrorCode PEPJDSetMinimalityIndex(PEP,PetscInt);
364: SLEPC_EXTERN PetscErrorCode PEPJDGetMinimalityIndex(PEP,PetscInt*);
365: SLEPC_EXTERN PetscErrorCode PEPJDSetProjection(PEP,PEPJDProjection);
366: SLEPC_EXTERN PetscErrorCode PEPJDGetProjection(PEP,PEPJDProjection*);

368: /*E
369:     PEPCISSExtraction - determines the extraction technique in the CISS solver

371:     Level: advanced

373: .seealso: PEPCISSSetExtraction(), PEPCISSGetExtraction()
374: E*/
375: typedef enum { PEP_CISS_EXTRACTION_RITZ,
376:                PEP_CISS_EXTRACTION_HANKEL,
377:                PEP_CISS_EXTRACTION_CAA    } PEPCISSExtraction;
378: SLEPC_EXTERN const char *PEPCISSExtractions[];

380: #if defined(PETSC_USE_COMPLEX) || defined(PETSC_CLANG_STATIC_ANALYZER)
381: SLEPC_EXTERN PetscErrorCode PEPCISSSetExtraction(PEP,PEPCISSExtraction);
382: SLEPC_EXTERN PetscErrorCode PEPCISSGetExtraction(PEP,PEPCISSExtraction*);
383: SLEPC_EXTERN PetscErrorCode PEPCISSSetSizes(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
384: SLEPC_EXTERN PetscErrorCode PEPCISSGetSizes(PEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
385: SLEPC_EXTERN PetscErrorCode PEPCISSSetThreshold(PEP,PetscReal,PetscReal);
386: SLEPC_EXTERN PetscErrorCode PEPCISSGetThreshold(PEP,PetscReal*,PetscReal*);
387: SLEPC_EXTERN PetscErrorCode PEPCISSSetRefinement(PEP,PetscInt,PetscInt);
388: SLEPC_EXTERN PetscErrorCode PEPCISSGetRefinement(PEP,PetscInt*,PetscInt*);
389: SLEPC_EXTERN PetscErrorCode PEPCISSGetKSPs(PEP,PetscInt*,KSP**);
390: #else
391: #define SlepcPEPCISSUnavailable(pep) do { \
393:     SETERRQ1(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"%s() not available with real scalars",PETSC_FUNCTION_NAME); \
394:     } while (0)
395: PETSC_STATIC_INLINE PetscErrorCode PEPCISSSetExtraction(PEP pep,PETSC_UNUSED PEPCISSExtraction ex) {SlepcPEPCISSUnavailable(pep);}
396: PETSC_STATIC_INLINE PetscErrorCode PEPCISSGetExtraction(PEP pep,PETSC_UNUSED PEPCISSExtraction *ex) {SlepcPEPCISSUnavailable(pep);}
397: PETSC_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);}
398: PETSC_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);}
399: PETSC_STATIC_INLINE PetscErrorCode PEPCISSSetThreshold(PEP pep,PETSC_UNUSED PetscReal delta,PETSC_UNUSED PetscReal spur) {SlepcPEPCISSUnavailable(pep);}
400: PETSC_STATIC_INLINE PetscErrorCode PEPCISSGetThreshold(PEP pep,PETSC_UNUSED PetscReal *delta,PETSC_UNUSED PetscReal *spur) {SlepcPEPCISSUnavailable(pep);}
401: PETSC_STATIC_INLINE PetscErrorCode PEPCISSSetRefinement(PEP pep,PETSC_UNUSED PetscInt inner,PETSC_UNUSED PetscInt blsize) {SlepcPEPCISSUnavailable(pep);}
402: PETSC_STATIC_INLINE PetscErrorCode PEPCISSGetRefinement(PEP pep,PETSC_UNUSED PetscInt *inner,PETSC_UNUSED PetscInt *blsize) {SlepcPEPCISSUnavailable(pep);}
403: PETSC_STATIC_INLINE PetscErrorCode PEPCISSGetKSPs(PEP pep,PETSC_UNUSED PetscInt *nsolve,PETSC_UNUSED KSP **ksp) {SlepcPEPCISSUnavailable(pep);}
404: #undef SlepcPEPCISSUnavailable
405: #endif
406: #endif