Actual source code: slepcnep.h
slepc-3.13.0 2020-03-31
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 nonlinear eigenvalue solvers
12: */
14: #if !defined(SLEPCNEP_H)
15: #define SLEPCNEP_H
16: #include <slepceps.h>
17: #include <slepcpep.h>
18: #include <slepcfn.h>
20: SLEPC_EXTERN PetscErrorCode NEPInitializePackage(void);
22: /*S
23: NEP - Abstract SLEPc object that manages all solvers for
24: nonlinear eigenvalue problems.
26: Level: beginner
28: .seealso: NEPCreate()
29: S*/
30: typedef struct _p_NEP* NEP;
32: /*J
33: NEPType - String with the name of a nonlinear eigensolver
35: Level: beginner
37: .seealso: NEPSetType(), NEP
38: J*/
39: typedef const char* NEPType;
40: #define NEPRII "rii"
41: #define NEPSLP "slp"
42: #define NEPNARNOLDI "narnoldi"
43: #define NEPCISS "ciss"
44: #define NEPINTERPOL "interpol"
45: #define NEPNLEIGS "nleigs"
47: /* Logging support */
48: SLEPC_EXTERN PetscClassId NEP_CLASSID;
50: /*E
51: NEPProblemType - Determines the type of the nonlinear eigenproblem
53: Level: intermediate
55: .seealso: NEPSetProblemType(), NEPGetProblemType()
56: E*/
57: typedef enum { NEP_GENERAL=1,
58: NEP_RATIONAL /* NEP defined in split form with all f_i rational */
59: } NEPProblemType;
61: /*E
62: NEPWhich - Determines which part of the spectrum is requested
64: Level: intermediate
66: .seealso: NEPSetWhichEigenpairs(), NEPGetWhichEigenpairs()
67: E*/
68: typedef enum { NEP_LARGEST_MAGNITUDE=1,
69: NEP_SMALLEST_MAGNITUDE,
70: NEP_LARGEST_REAL,
71: NEP_SMALLEST_REAL,
72: NEP_LARGEST_IMAGINARY,
73: NEP_SMALLEST_IMAGINARY,
74: NEP_TARGET_MAGNITUDE,
75: NEP_TARGET_REAL,
76: NEP_TARGET_IMAGINARY,
77: NEP_ALL,
78: NEP_WHICH_USER } NEPWhich;
80: /*E
81: NEPErrorType - The error type used to assess accuracy of computed solutions
83: Level: intermediate
85: .seealso: NEPComputeError()
86: E*/
87: typedef enum { NEP_ERROR_ABSOLUTE,
88: NEP_ERROR_RELATIVE,
89: NEP_ERROR_BACKWARD } NEPErrorType;
90: SLEPC_EXTERN const char *NEPErrorTypes[];
92: /*E
93: NEPRefine - The refinement type
95: Level: intermediate
97: .seealso: NEPSetRefine()
98: E*/
99: typedef enum { NEP_REFINE_NONE,
100: NEP_REFINE_SIMPLE,
101: NEP_REFINE_MULTIPLE } NEPRefine;
102: SLEPC_EXTERN const char *NEPRefineTypes[];
104: /*E
105: NEPRefineScheme - The scheme used for solving linear systems during iterative refinement
107: Level: intermediate
109: .seealso: NEPSetRefine()
110: E*/
111: typedef enum { NEP_REFINE_SCHEME_SCHUR=1,
112: NEP_REFINE_SCHEME_MBE,
113: NEP_REFINE_SCHEME_EXPLICIT } NEPRefineScheme;
114: SLEPC_EXTERN const char *NEPRefineSchemes[];
116: /*E
117: NEPConv - Determines the convergence test
119: Level: intermediate
121: .seealso: NEPSetConvergenceTest(), NEPSetConvergenceTestFunction()
122: E*/
123: typedef enum { NEP_CONV_ABS,
124: NEP_CONV_REL,
125: NEP_CONV_NORM,
126: NEP_CONV_USER } NEPConv;
128: /*E
129: NEPStop - Determines the stopping test
131: Level: advanced
133: .seealso: NEPSetStoppingTest(), NEPSetStoppingTestFunction()
134: E*/
135: typedef enum { NEP_STOP_BASIC,
136: NEP_STOP_USER } NEPStop;
138: /*E
139: NEPConvergedReason - Reason a nonlinear eigensolver was said to
140: have converged or diverged
142: Level: intermediate
144: .seealso: NEPSolve(), NEPGetConvergedReason(), NEPSetTolerances()
145: E*/
146: typedef enum {/* converged */
147: NEP_CONVERGED_TOL = 1,
148: NEP_CONVERGED_USER = 2,
149: /* diverged */
150: NEP_DIVERGED_ITS = -1,
151: NEP_DIVERGED_BREAKDOWN = -2,
152: /* unused = -3 */
153: NEP_DIVERGED_LINEAR_SOLVE = -4,
154: NEP_DIVERGED_SUBSPACE_EXHAUSTED = -5,
155: NEP_CONVERGED_ITERATING = 0} NEPConvergedReason;
156: SLEPC_EXTERN const char *const*NEPConvergedReasons;
158: SLEPC_EXTERN PetscErrorCode NEPCreate(MPI_Comm,NEP*);
159: SLEPC_EXTERN PetscErrorCode NEPDestroy(NEP*);
160: SLEPC_EXTERN PetscErrorCode NEPReset(NEP);
161: SLEPC_EXTERN PetscErrorCode NEPSetType(NEP,NEPType);
162: SLEPC_EXTERN PetscErrorCode NEPGetType(NEP,NEPType*);
163: SLEPC_EXTERN PetscErrorCode NEPSetProblemType(NEP,NEPProblemType);
164: SLEPC_EXTERN PetscErrorCode NEPGetProblemType(NEP,NEPProblemType*);
165: SLEPC_EXTERN PetscErrorCode NEPSetTarget(NEP,PetscScalar);
166: SLEPC_EXTERN PetscErrorCode NEPGetTarget(NEP,PetscScalar*);
167: SLEPC_EXTERN PetscErrorCode NEPSetFromOptions(NEP);
168: SLEPC_EXTERN PetscErrorCode NEPSetUp(NEP);
169: SLEPC_EXTERN PetscErrorCode NEPSolve(NEP);
170: SLEPC_EXTERN PetscErrorCode NEPView(NEP,PetscViewer);
171: SLEPC_EXTERN PetscErrorCode NEPViewFromOptions(NEP,PetscObject,const char[]);
172: SLEPC_EXTERN PetscErrorCode NEPErrorView(NEP,NEPErrorType,PetscViewer);
173: SLEPC_EXTERN PetscErrorCode NEPErrorViewFromOptions(NEP);
174: SLEPC_EXTERN PetscErrorCode NEPReasonView(NEP,PetscViewer);
175: SLEPC_EXTERN PetscErrorCode NEPReasonViewFromOptions(NEP);
176: SLEPC_EXTERN PetscErrorCode NEPValuesView(NEP,PetscViewer);
177: SLEPC_EXTERN PetscErrorCode NEPValuesViewFromOptions(NEP);
178: SLEPC_EXTERN PetscErrorCode NEPVectorsView(NEP,PetscViewer);
179: SLEPC_EXTERN PetscErrorCode NEPVectorsViewFromOptions(NEP);
181: SLEPC_EXTERN PetscErrorCode NEPSetFunction(NEP,Mat,Mat,PetscErrorCode (*)(NEP,PetscScalar,Mat,Mat,void*),void*);
182: SLEPC_EXTERN PetscErrorCode NEPGetFunction(NEP,Mat*,Mat*,PetscErrorCode (**)(NEP,PetscScalar,Mat,Mat,void*),void**);
183: SLEPC_EXTERN PetscErrorCode NEPSetJacobian(NEP,Mat,PetscErrorCode (*)(NEP,PetscScalar,Mat,void*),void*);
184: SLEPC_EXTERN PetscErrorCode NEPGetJacobian(NEP,Mat*,PetscErrorCode (**)(NEP,PetscScalar,Mat,void*),void**);
185: PETSC_DEPRECATED_FUNCTION("Use NEPSetFunction() and NEPSetJacobian") PETSC_STATIC_INLINE PetscErrorCode NEPSetDerivatives(NEP nep,Mat A,PetscErrorCode (*fun)(NEP,PetscScalar,PetscInt,Mat,void*),void *ctx) {SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented in this version");}
186: PETSC_DEPRECATED_FUNCTION("Use NEPGetFunction() and NEPGetJacobian") PETSC_STATIC_INLINE PetscErrorCode NEPGetDerivatives(NEP nep,Mat *A,PetscErrorCode (**fun)(NEP,PetscScalar,PetscInt,Mat,void*),void **ctx) {SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented in this version");}
187: SLEPC_EXTERN PetscErrorCode NEPSetSplitOperator(NEP,PetscInt,Mat[],FN[],MatStructure);
188: SLEPC_EXTERN PetscErrorCode NEPGetSplitOperatorTerm(NEP,PetscInt,Mat*,FN*);
189: SLEPC_EXTERN PetscErrorCode NEPGetSplitOperatorInfo(NEP,PetscInt*,MatStructure*);
191: SLEPC_EXTERN PetscErrorCode NEPSetBV(NEP,BV);
192: SLEPC_EXTERN PetscErrorCode NEPGetBV(NEP,BV*);
193: SLEPC_EXTERN PetscErrorCode NEPSetRG(NEP,RG);
194: SLEPC_EXTERN PetscErrorCode NEPGetRG(NEP,RG*);
195: SLEPC_EXTERN PetscErrorCode NEPSetDS(NEP,DS);
196: SLEPC_EXTERN PetscErrorCode NEPGetDS(NEP,DS*);
197: SLEPC_EXTERN PetscErrorCode NEPRefineGetKSP(NEP,KSP*);
198: SLEPC_EXTERN PetscErrorCode NEPSetTolerances(NEP,PetscReal,PetscInt);
199: SLEPC_EXTERN PetscErrorCode NEPGetTolerances(NEP,PetscReal*,PetscInt*);
200: SLEPC_EXTERN PetscErrorCode NEPSetConvergenceTestFunction(NEP,PetscErrorCode (*)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void*,PetscErrorCode (*)(void*));
201: SLEPC_EXTERN PetscErrorCode NEPSetConvergenceTest(NEP,NEPConv);
202: SLEPC_EXTERN PetscErrorCode NEPGetConvergenceTest(NEP,NEPConv*);
203: SLEPC_EXTERN PetscErrorCode NEPConvergedAbsolute(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
204: SLEPC_EXTERN PetscErrorCode NEPConvergedRelative(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
205: SLEPC_EXTERN PetscErrorCode NEPConvergedNorm(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
206: SLEPC_EXTERN PetscErrorCode NEPSetStoppingTestFunction(NEP,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
207: SLEPC_EXTERN PetscErrorCode NEPSetStoppingTest(NEP,NEPStop);
208: SLEPC_EXTERN PetscErrorCode NEPGetStoppingTest(NEP,NEPStop*);
209: SLEPC_EXTERN PetscErrorCode NEPStoppingBasic(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);
210: SLEPC_EXTERN PetscErrorCode NEPSetDimensions(NEP,PetscInt,PetscInt,PetscInt);
211: SLEPC_EXTERN PetscErrorCode NEPGetDimensions(NEP,PetscInt*,PetscInt*,PetscInt*);
212: SLEPC_EXTERN PetscErrorCode NEPSetRefine(NEP,NEPRefine,PetscInt,PetscReal,PetscInt,NEPRefineScheme);
213: SLEPC_EXTERN PetscErrorCode NEPGetRefine(NEP,NEPRefine*,PetscInt*,PetscReal*,PetscInt*,NEPRefineScheme*);
215: SLEPC_EXTERN PetscErrorCode NEPGetConverged(NEP,PetscInt*);
216: SLEPC_EXTERN PetscErrorCode NEPGetEigenpair(NEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
217: SLEPC_EXTERN PetscErrorCode NEPGetLeftEigenvector(NEP,PetscInt,Vec,Vec);
219: SLEPC_EXTERN PetscErrorCode NEPComputeError(NEP,PetscInt,NEPErrorType,PetscReal*);
220: PETSC_DEPRECATED_FUNCTION("Use NEPComputeError()") PETSC_STATIC_INLINE PetscErrorCode NEPComputeRelativeError(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_RELATIVE,r);}
221: PETSC_DEPRECATED_FUNCTION("Use NEPComputeError() with NEP_ERROR_ABSOLUTE") PETSC_STATIC_INLINE PetscErrorCode NEPComputeResidualNorm(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_ABSOLUTE,r);}
222: SLEPC_EXTERN PetscErrorCode NEPGetErrorEstimate(NEP,PetscInt,PetscReal*);
224: SLEPC_EXTERN PetscErrorCode NEPComputeFunction(NEP,PetscScalar,Mat,Mat);
225: SLEPC_EXTERN PetscErrorCode NEPComputeJacobian(NEP,PetscScalar,Mat);
226: SLEPC_EXTERN PetscErrorCode NEPApplyFunction(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
227: SLEPC_EXTERN PetscErrorCode NEPApplyAdjoint(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
228: SLEPC_EXTERN PetscErrorCode NEPApplyJacobian(NEP,PetscScalar,Vec,Vec,Vec,Mat);
229: SLEPC_EXTERN PetscErrorCode NEPProjectOperator(NEP,PetscInt,PetscInt);
231: SLEPC_EXTERN PetscErrorCode NEPMonitor(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
232: SLEPC_EXTERN PetscErrorCode NEPMonitorSet(NEP,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
233: SLEPC_EXTERN PetscErrorCode NEPMonitorSetFromOptions(NEP,const char*,const char*,const char*,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool);
234: SLEPC_EXTERN PetscErrorCode NEPConvMonitorSetFromOptions(NEP,const char*,const char*,const char*,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor));
235: SLEPC_EXTERN PetscErrorCode NEPMonitorCancel(NEP);
236: SLEPC_EXTERN PetscErrorCode NEPGetMonitorContext(NEP,void **);
237: SLEPC_EXTERN PetscErrorCode NEPGetIterationNumber(NEP,PetscInt*);
239: SLEPC_EXTERN PetscErrorCode NEPSetInitialSpace(NEP,PetscInt,Vec[]);
240: SLEPC_EXTERN PetscErrorCode NEPSetWhichEigenpairs(NEP,NEPWhich);
241: SLEPC_EXTERN PetscErrorCode NEPGetWhichEigenpairs(NEP,NEPWhich*);
242: SLEPC_EXTERN PetscErrorCode NEPSetTwoSided(NEP,PetscBool);
243: SLEPC_EXTERN PetscErrorCode NEPGetTwoSided(NEP,PetscBool*);
245: SLEPC_EXTERN PetscErrorCode NEPApplyResolvent(NEP,RG,PetscScalar,Vec,Vec);
247: SLEPC_EXTERN PetscErrorCode NEPSetEigenvalueComparison(NEP,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);
249: SLEPC_EXTERN PetscErrorCode NEPMonitorAll(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
250: SLEPC_EXTERN PetscErrorCode NEPMonitorFirst(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
251: SLEPC_EXTERN PetscErrorCode NEPMonitorConverged(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor);
252: SLEPC_EXTERN PetscErrorCode NEPMonitorLGCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
253: SLEPC_EXTERN PetscErrorCode NEPMonitorLG(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
254: SLEPC_EXTERN PetscErrorCode NEPMonitorLGAll(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
256: SLEPC_EXTERN PetscErrorCode NEPSetTrackAll(NEP,PetscBool);
257: SLEPC_EXTERN PetscErrorCode NEPGetTrackAll(NEP,PetscBool*);
259: SLEPC_EXTERN PetscErrorCode NEPSetOptionsPrefix(NEP,const char*);
260: SLEPC_EXTERN PetscErrorCode NEPAppendOptionsPrefix(NEP,const char*);
261: SLEPC_EXTERN PetscErrorCode NEPGetOptionsPrefix(NEP,const char*[]);
263: SLEPC_EXTERN PetscErrorCode NEPGetConvergedReason(NEP,NEPConvergedReason*);
265: SLEPC_EXTERN PetscFunctionList NEPList;
266: SLEPC_EXTERN PetscErrorCode NEPRegister(const char[],PetscErrorCode(*)(NEP));
268: SLEPC_EXTERN PetscErrorCode NEPSetWorkVecs(NEP,PetscInt);
269: SLEPC_EXTERN PetscErrorCode NEPAllocateSolution(NEP,PetscInt);
271: /* --------- options specific to particular eigensolvers -------- */
273: SLEPC_EXTERN PetscErrorCode NEPRIISetMaximumIterations(NEP,PetscInt);
274: SLEPC_EXTERN PetscErrorCode NEPRIIGetMaximumIterations(NEP,PetscInt*);
275: SLEPC_EXTERN PetscErrorCode NEPRIISetLagPreconditioner(NEP,PetscInt);
276: SLEPC_EXTERN PetscErrorCode NEPRIIGetLagPreconditioner(NEP,PetscInt*);
277: SLEPC_EXTERN PetscErrorCode NEPRIISetConstCorrectionTol(NEP,PetscBool);
278: SLEPC_EXTERN PetscErrorCode NEPRIIGetConstCorrectionTol(NEP,PetscBool*);
279: SLEPC_EXTERN PetscErrorCode NEPRIISetHermitian(NEP,PetscBool);
280: SLEPC_EXTERN PetscErrorCode NEPRIIGetHermitian(NEP,PetscBool*);
281: SLEPC_EXTERN PetscErrorCode NEPRIISetDeflationThreshold(NEP,PetscReal);
282: SLEPC_EXTERN PetscErrorCode NEPRIIGetDeflationThreshold(NEP,PetscReal*);
283: SLEPC_EXTERN PetscErrorCode NEPRIISetKSP(NEP,KSP);
284: SLEPC_EXTERN PetscErrorCode NEPRIIGetKSP(NEP,KSP*);
286: SLEPC_EXTERN PetscErrorCode NEPSLPSetDeflationThreshold(NEP,PetscReal);
287: SLEPC_EXTERN PetscErrorCode NEPSLPGetDeflationThreshold(NEP,PetscReal*);
288: SLEPC_EXTERN PetscErrorCode NEPSLPSetEPS(NEP,EPS);
289: SLEPC_EXTERN PetscErrorCode NEPSLPGetEPS(NEP,EPS*);
290: SLEPC_EXTERN PetscErrorCode NEPSLPSetEPSLeft(NEP,EPS);
291: SLEPC_EXTERN PetscErrorCode NEPSLPGetEPSLeft(NEP,EPS*);
292: SLEPC_EXTERN PetscErrorCode NEPSLPSetKSP(NEP,KSP);
293: SLEPC_EXTERN PetscErrorCode NEPSLPGetKSP(NEP,KSP*);
295: SLEPC_EXTERN PetscErrorCode NEPNArnoldiSetKSP(NEP,KSP);
296: SLEPC_EXTERN PetscErrorCode NEPNArnoldiGetKSP(NEP,KSP*);
297: SLEPC_EXTERN PetscErrorCode NEPNArnoldiSetLagPreconditioner(NEP,PetscInt);
298: SLEPC_EXTERN PetscErrorCode NEPNArnoldiGetLagPreconditioner(NEP,PetscInt*);
300: #if defined(PETSC_USE_COMPLEX)
301: SLEPC_EXTERN PetscErrorCode NEPCISSSetSizes(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
302: SLEPC_EXTERN PetscErrorCode NEPCISSGetSizes(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
303: SLEPC_EXTERN PetscErrorCode NEPCISSSetThreshold(NEP,PetscReal,PetscReal);
304: SLEPC_EXTERN PetscErrorCode NEPCISSGetThreshold(NEP,PetscReal*,PetscReal*);
305: SLEPC_EXTERN PetscErrorCode NEPCISSSetRefinement(NEP,PetscInt,PetscInt);
306: SLEPC_EXTERN PetscErrorCode NEPCISSGetRefinement(NEP,PetscInt*,PetscInt*);
307: SLEPC_EXTERN PetscErrorCode NEPCISSGetKSPs(NEP,PetscInt*,KSP**);
308: #else
309: #define SlepcNEPCISSUnavailable(nep) do { \
311: SETERRQ1(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"%s() not available with real scalars",PETSC_FUNCTION_NAME); \
312: } while (0)
313: PETSC_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);}
314: PETSC_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);}
315: PETSC_STATIC_INLINE PetscErrorCode NEPCISSSetThreshold(NEP nep,PETSC_UNUSED PetscReal delta,PETSC_UNUSED PetscReal spur) {SlepcNEPCISSUnavailable(nep);}
316: PETSC_STATIC_INLINE PetscErrorCode NEPCISSGetThreshold(NEP nep,PETSC_UNUSED PetscReal *delta,PETSC_UNUSED PetscReal *spur) {SlepcNEPCISSUnavailable(nep);}
317: PETSC_STATIC_INLINE PetscErrorCode NEPCISSSetRefinement(NEP nep,PETSC_UNUSED PetscInt inner,PETSC_UNUSED PetscInt blsize) {SlepcNEPCISSUnavailable(nep);}
318: PETSC_STATIC_INLINE PetscErrorCode NEPCISSGetRefinement(NEP nep,PETSC_UNUSED PetscInt *inner,PETSC_UNUSED PetscInt *blsize) {SlepcNEPCISSUnavailable(nep);}
319: PETSC_STATIC_INLINE PetscErrorCode NEPCISSGetKSPs(NEP nep,PETSC_UNUSED PetscInt *nsolve,PETSC_UNUSED KSP **ksp) {SlepcNEPCISSUnavailable(nep);}
320: #undef SlepcNEPCISSUnavailable
321: #endif
323: SLEPC_EXTERN PetscErrorCode NEPInterpolSetPEP(NEP,PEP);
324: SLEPC_EXTERN PetscErrorCode NEPInterpolGetPEP(NEP,PEP*);
325: SLEPC_EXTERN PetscErrorCode NEPInterpolSetInterpolation(NEP,PetscReal,PetscInt);
326: SLEPC_EXTERN PetscErrorCode NEPInterpolGetInterpolation(NEP,PetscReal*,PetscInt*);
328: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetSingularitiesFunction(NEP,PetscErrorCode (*)(NEP,PetscInt*,PetscScalar*,void*),void*);
329: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetSingularitiesFunction(NEP,PetscErrorCode (**)(NEP,PetscInt*,PetscScalar*,void*),void **);
330: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetRestart(NEP,PetscReal);
331: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetRestart(NEP,PetscReal*);
332: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetLocking(NEP,PetscBool);
333: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetLocking(NEP,PetscBool*);
334: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetInterpolation(NEP,PetscReal,PetscInt);
335: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetInterpolation(NEP,PetscReal*,PetscInt*);
336: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetRKShifts(NEP,PetscInt,PetscScalar[]);
337: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetRKShifts(NEP,PetscInt*,PetscScalar*[]);
338: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetKSPs(NEP,PetscInt*,KSP**);
339: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetFullBasis(NEP,PetscBool);
340: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetFullBasis(NEP,PetscBool*);
341: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetEPS(NEP,EPS);
342: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetEPS(NEP,EPS*);
344: #endif