Actual source code: slepcpep.h

slepc-3.14.0 2020-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 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"

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

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

 50:     Level: intermediate

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

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

 63:     Level: intermediate

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

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

 83:     Level: intermediate

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

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

 98:     Level: intermediate

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

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

111:     Level: intermediate

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

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

123:     Level: intermediate

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

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

135:     Level: intermediate

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

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

148:     Level: intermediate

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

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

160:     Level: intermediate

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

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

172:     Level: advanced

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

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

183:     Level: intermediate

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

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

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

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

262: SLEPC_EXTERN PetscErrorCode PEPGetConverged(PEP,PetscInt*);
263: SLEPC_EXTERN PetscErrorCode PEPGetEigenpair(PEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
264: SLEPC_EXTERN PetscErrorCode PEPComputeError(PEP,PetscInt,PEPErrorType,PetscReal*);
265: PETSC_DEPRECATED_FUNCTION("Use PEPComputeError()") PETSC_STATIC_INLINE PetscErrorCode PEPComputeRelativeError(PEP pep,PetscInt i,PetscReal *r) {return PEPComputeError(pep,i,PEP_ERROR_BACKWARD,r);}
266: 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);}
267: SLEPC_EXTERN PetscErrorCode PEPGetErrorEstimate(PEP,PetscInt,PetscReal*);

269: SLEPC_EXTERN PetscErrorCode PEPMonitor(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
270: SLEPC_EXTERN PetscErrorCode PEPMonitorSet(PEP,PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
271: SLEPC_EXTERN PetscErrorCode PEPMonitorSetFromOptions(PEP,const char*,const char*,const char*,PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool);
272: SLEPC_EXTERN PetscErrorCode PEPConvMonitorSetFromOptions(PEP,const char*,const char*,const char*,PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor));
273: SLEPC_EXTERN PetscErrorCode PEPMonitorCancel(PEP);
274: SLEPC_EXTERN PetscErrorCode PEPGetMonitorContext(PEP,void **);
275: SLEPC_EXTERN PetscErrorCode PEPGetIterationNumber(PEP,PetscInt*);

277: SLEPC_EXTERN PetscErrorCode PEPSetInitialSpace(PEP,PetscInt,Vec[]);
278: SLEPC_EXTERN PetscErrorCode PEPSetWhichEigenpairs(PEP,PEPWhich);
279: SLEPC_EXTERN PetscErrorCode PEPGetWhichEigenpairs(PEP,PEPWhich*);
280: SLEPC_EXTERN PetscErrorCode PEPSetEigenvalueComparison(PEP,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);

282: SLEPC_EXTERN PetscErrorCode PEPMonitorAll(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
283: SLEPC_EXTERN PetscErrorCode PEPMonitorFirst(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
284: SLEPC_EXTERN PetscErrorCode PEPMonitorConverged(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor);
285: SLEPC_EXTERN PetscErrorCode PEPMonitorLGCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
286: SLEPC_EXTERN PetscErrorCode PEPMonitorLG(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
287: SLEPC_EXTERN PetscErrorCode PEPMonitorLGAll(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);

289: SLEPC_EXTERN PetscErrorCode PEPSetTrackAll(PEP,PetscBool);
290: SLEPC_EXTERN PetscErrorCode PEPGetTrackAll(PEP,PetscBool*);

292: SLEPC_EXTERN PetscErrorCode PEPSetOptionsPrefix(PEP,const char*);
293: SLEPC_EXTERN PetscErrorCode PEPAppendOptionsPrefix(PEP,const char*);
294: SLEPC_EXTERN PetscErrorCode PEPGetOptionsPrefix(PEP,const char*[]);

296: SLEPC_EXTERN PetscFunctionList PEPList;
297: SLEPC_EXTERN PetscErrorCode PEPRegister(const char[],PetscErrorCode(*)(PEP));

299: SLEPC_EXTERN PetscErrorCode PEPSetWorkVecs(PEP,PetscInt);
300: SLEPC_EXTERN PetscErrorCode PEPAllocateSolution(PEP,PetscInt);

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

304: SLEPC_EXTERN PetscErrorCode PEPLinearSetLinearization(PEP,PetscReal,PetscReal);
305: SLEPC_EXTERN PetscErrorCode PEPLinearGetLinearization(PEP,PetscReal*,PetscReal*);
306: SLEPC_EXTERN PetscErrorCode PEPLinearSetExplicitMatrix(PEP,PetscBool);
307: SLEPC_EXTERN PetscErrorCode PEPLinearGetExplicitMatrix(PEP,PetscBool*);
308: SLEPC_EXTERN PetscErrorCode PEPLinearSetEPS(PEP,EPS);
309: SLEPC_EXTERN PetscErrorCode PEPLinearGetEPS(PEP,EPS*);
310: 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);}
311: PETSC_DEPRECATED_FUNCTION("Use PEPLinearGetLinearization()") PETSC_STATIC_INLINE PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform) {(void)pep; if (cform) *cform=1; return 0;}

313: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetRestart(PEP,PetscReal);
314: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetRestart(PEP,PetscReal*);
315: SLEPC_EXTERN PetscErrorCode PEPQArnoldiSetLocking(PEP,PetscBool);
316: SLEPC_EXTERN PetscErrorCode PEPQArnoldiGetLocking(PEP,PetscBool*);

318: SLEPC_EXTERN PetscErrorCode PEPTOARSetRestart(PEP,PetscReal);
319: SLEPC_EXTERN PetscErrorCode PEPTOARGetRestart(PEP,PetscReal*);
320: SLEPC_EXTERN PetscErrorCode PEPTOARSetLocking(PEP,PetscBool);
321: SLEPC_EXTERN PetscErrorCode PEPTOARGetLocking(PEP,PetscBool*);

323: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLinearization(PEP,PetscReal,PetscReal);
324: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLinearization(PEP,PetscReal*,PetscReal*);
325: SLEPC_EXTERN PetscErrorCode PEPSTOARSetLocking(PEP,PetscBool);
326: SLEPC_EXTERN PetscErrorCode PEPSTOARGetLocking(PEP,PetscBool*);
327: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDetectZeros(PEP,PetscBool);
328: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDetectZeros(PEP,PetscBool*);
329: SLEPC_EXTERN PetscErrorCode PEPSTOARGetInertias(PEP,PetscInt*,PetscReal**,PetscInt**);
330: SLEPC_EXTERN PetscErrorCode PEPSTOARSetDimensions(PEP,PetscInt,PetscInt,PetscInt);
331: SLEPC_EXTERN PetscErrorCode PEPSTOARGetDimensions(PEP,PetscInt*,PetscInt*,PetscInt*);
332: SLEPC_EXTERN PetscErrorCode PEPSTOARSetCheckEigenvalueType(PEP,PetscBool);
333: SLEPC_EXTERN PetscErrorCode PEPSTOARGetCheckEigenvalueType(PEP,PetscBool*);
334: SLEPC_EXTERN PetscErrorCode PEPCheckDefiniteQEP(PEP,PetscReal*,PetscReal*,PetscInt*,PetscInt*);

336: /*E
337:     PEPJDProjection - The scheme used for solving linear systems during iterative refinement

339:     Level: intermediate

341: .seealso: PEPJDSetProjection()
342: E*/
343: typedef enum { PEP_JD_PROJECTION_HARMONIC,
344:                PEP_JD_PROJECTION_ORTHOGONAL } PEPJDProjection;
345: SLEPC_EXTERN const char *PEPJDProjectionTypes[];

347: SLEPC_EXTERN PetscErrorCode PEPJDSetRestart(PEP,PetscReal);
348: SLEPC_EXTERN PetscErrorCode PEPJDGetRestart(PEP,PetscReal*);
349: SLEPC_EXTERN PetscErrorCode PEPJDSetFix(PEP,PetscReal);
350: SLEPC_EXTERN PetscErrorCode PEPJDGetFix(PEP,PetscReal*);
351: SLEPC_EXTERN PetscErrorCode PEPJDSetReusePreconditioner(PEP,PetscBool);
352: SLEPC_EXTERN PetscErrorCode PEPJDGetReusePreconditioner(PEP,PetscBool*);
353: SLEPC_EXTERN PetscErrorCode PEPJDSetMinimalityIndex(PEP,PetscInt);
354: SLEPC_EXTERN PetscErrorCode PEPJDGetMinimalityIndex(PEP,PetscInt*);
355: SLEPC_EXTERN PetscErrorCode PEPJDSetProjection(PEP,PEPJDProjection);
356: SLEPC_EXTERN PetscErrorCode PEPJDGetProjection(PEP,PEPJDProjection*);

358: #endif