Actual source code: slepcnep.h
slepc-main 2024-12-17
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 nonlinear eigenvalue solvers
12: */
14: #pragma once
16: #include <slepceps.h>
17: #include <slepcpep.h>
18: #include <slepcfn.h>
20: /* SUBMANSEC = NEP */
22: SLEPC_EXTERN PetscErrorCode NEPInitializePackage(void);
23: SLEPC_EXTERN PetscErrorCode NEPFinalizePackage(void);
25: /*S
26: NEP - Abstract SLEPc object that manages all solvers for
27: nonlinear eigenvalue problems.
29: Level: beginner
31: .seealso: NEPCreate()
32: S*/
33: typedef struct _p_NEP* NEP;
35: /*J
36: NEPType - String with the name of a nonlinear eigensolver
38: Level: beginner
40: .seealso: NEPSetType(), NEP
41: J*/
42: typedef const char* NEPType;
43: #define NEPRII "rii"
44: #define NEPSLP "slp"
45: #define NEPNARNOLDI "narnoldi"
46: #define NEPCISS "ciss"
47: #define NEPINTERPOL "interpol"
48: #define NEPNLEIGS "nleigs"
50: /* Logging support */
51: SLEPC_EXTERN PetscClassId NEP_CLASSID;
53: /*E
54: NEPProblemType - Determines the type of the nonlinear eigenproblem
56: Level: intermediate
58: .seealso: NEPSetProblemType(), NEPGetProblemType()
59: E*/
60: typedef enum { NEP_GENERAL=1,
61: NEP_RATIONAL /* NEP defined in split form with all f_i rational */
62: } NEPProblemType;
64: /*E
65: NEPWhich - Determines which part of the spectrum is requested
67: Level: intermediate
69: .seealso: NEPSetWhichEigenpairs(), NEPGetWhichEigenpairs()
70: E*/
71: typedef enum { NEP_LARGEST_MAGNITUDE=1,
72: NEP_SMALLEST_MAGNITUDE,
73: NEP_LARGEST_REAL,
74: NEP_SMALLEST_REAL,
75: NEP_LARGEST_IMAGINARY,
76: NEP_SMALLEST_IMAGINARY,
77: NEP_TARGET_MAGNITUDE,
78: NEP_TARGET_REAL,
79: NEP_TARGET_IMAGINARY,
80: NEP_ALL,
81: NEP_WHICH_USER } NEPWhich;
83: /*E
84: NEPErrorType - The error type used to assess accuracy of computed solutions
86: Level: intermediate
88: .seealso: NEPComputeError()
89: E*/
90: typedef enum { NEP_ERROR_ABSOLUTE,
91: NEP_ERROR_RELATIVE,
92: NEP_ERROR_BACKWARD } NEPErrorType;
93: SLEPC_EXTERN const char *NEPErrorTypes[];
95: /*E
96: NEPRefine - The refinement type
98: Level: intermediate
100: .seealso: NEPSetRefine()
101: E*/
102: typedef enum { NEP_REFINE_NONE,
103: NEP_REFINE_SIMPLE,
104: NEP_REFINE_MULTIPLE } NEPRefine;
105: SLEPC_EXTERN const char *NEPRefineTypes[];
107: /*E
108: NEPRefineScheme - The scheme used for solving linear systems during iterative refinement
110: Level: intermediate
112: .seealso: NEPSetRefine()
113: E*/
114: typedef enum { NEP_REFINE_SCHEME_SCHUR=1,
115: NEP_REFINE_SCHEME_MBE,
116: NEP_REFINE_SCHEME_EXPLICIT } NEPRefineScheme;
117: SLEPC_EXTERN const char *NEPRefineSchemes[];
119: /*E
120: NEPConv - Determines the convergence test
122: Level: intermediate
124: .seealso: NEPSetConvergenceTest(), NEPSetConvergenceTestFunction()
125: E*/
126: typedef enum { NEP_CONV_ABS,
127: NEP_CONV_REL,
128: NEP_CONV_NORM,
129: NEP_CONV_USER } NEPConv;
131: /*E
132: NEPStop - Determines the stopping test
134: Level: advanced
136: .seealso: NEPSetStoppingTest(), NEPSetStoppingTestFunction()
137: E*/
138: typedef enum { NEP_STOP_BASIC,
139: NEP_STOP_USER } NEPStop;
141: /*E
142: NEPConvergedReason - Reason a nonlinear eigensolver was said to
143: have converged or diverged
145: Level: intermediate
147: .seealso: NEPSolve(), NEPGetConvergedReason(), NEPSetTolerances()
148: E*/
149: typedef enum {/* converged */
150: NEP_CONVERGED_TOL = 1,
151: NEP_CONVERGED_USER = 2,
152: /* diverged */
153: NEP_DIVERGED_ITS = -1,
154: NEP_DIVERGED_BREAKDOWN = -2,
155: /* unused = -3 */
156: NEP_DIVERGED_LINEAR_SOLVE = -4,
157: NEP_DIVERGED_SUBSPACE_EXHAUSTED = -5,
158: NEP_CONVERGED_ITERATING = 0} NEPConvergedReason;
159: SLEPC_EXTERN const char *const*NEPConvergedReasons;
161: SLEPC_EXTERN PetscErrorCode NEPCreate(MPI_Comm,NEP*);
162: SLEPC_EXTERN PetscErrorCode NEPDestroy(NEP*);
163: SLEPC_EXTERN PetscErrorCode NEPReset(NEP);
164: SLEPC_EXTERN PetscErrorCode NEPSetType(NEP,NEPType);
165: SLEPC_EXTERN PetscErrorCode NEPGetType(NEP,NEPType*);
166: SLEPC_EXTERN PetscErrorCode NEPSetProblemType(NEP,NEPProblemType);
167: SLEPC_EXTERN PetscErrorCode NEPGetProblemType(NEP,NEPProblemType*);
168: SLEPC_EXTERN PetscErrorCode NEPSetTarget(NEP,PetscScalar);
169: SLEPC_EXTERN PetscErrorCode NEPGetTarget(NEP,PetscScalar*);
170: SLEPC_EXTERN PetscErrorCode NEPSetFromOptions(NEP);
171: SLEPC_EXTERN PetscErrorCode NEPSetDSType(NEP);
172: SLEPC_EXTERN PetscErrorCode NEPSetUp(NEP);
173: SLEPC_EXTERN PetscErrorCode NEPSolve(NEP);
174: SLEPC_EXTERN PetscErrorCode NEPView(NEP,PetscViewer);
175: SLEPC_EXTERN PetscErrorCode NEPViewFromOptions(NEP,PetscObject,const char[]);
176: SLEPC_EXTERN PetscErrorCode NEPErrorView(NEP,NEPErrorType,PetscViewer);
177: SLEPC_EXTERN PetscErrorCode NEPErrorViewFromOptions(NEP);
178: SLEPC_EXTERN PetscErrorCode NEPConvergedReasonView(NEP,PetscViewer);
179: SLEPC_EXTERN PetscErrorCode NEPConvergedReasonViewFromOptions(NEP);
180: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "NEPConvergedReasonView()", ) static inline PetscErrorCode NEPReasonView(NEP nep,PetscViewer v) {return NEPConvergedReasonView(nep,v);}
181: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "NEPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode NEPReasonViewFromOptions(NEP nep) {return NEPConvergedReasonViewFromOptions(nep);}
182: SLEPC_EXTERN PetscErrorCode NEPValuesView(NEP,PetscViewer);
183: SLEPC_EXTERN PetscErrorCode NEPValuesViewFromOptions(NEP);
184: SLEPC_EXTERN PetscErrorCode NEPVectorsView(NEP,PetscViewer);
185: SLEPC_EXTERN PetscErrorCode NEPVectorsViewFromOptions(NEP);
187: /*S
188: NEPFunctionFn - A prototype of a NEP function evaluation function that would be passed to NEPSetFunction()
190: Calling Sequence:
191: + nep - eigensolver context obtained from NEPCreate()
192: . lambda - the scalar argument where T(.) must be evaluated
193: . T - matrix that will contain T(lambda)
194: . P - [optional] different matrix to build the preconditioner
195: - ctx - [optional] user-defined context for private data for the
196: function evaluation routine (may be NULL)
198: Level: beginner
200: .seealso: NEPSetFunction()
201: S*/
202: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(NEPFunctionFn)(NEP nep,PetscScalar lambda,Mat T,Mat P,void *ctx);
204: /*S
205: NEPJacobianFn - A prototype of a NEP Jacobian evaluation function that would be passed to NEPSetJacobian()
207: Calling Sequence:
208: + nep - eigensolver context obtained from NEPCreate()
209: . lambda - the scalar argument where T'(.) must be evaluated
210: . J - matrix that will contain T'(lambda)
211: - ctx - [optional] user-defined context for private data for the
212: Jacobian evaluation routine (may be NULL)
214: Level: beginner
216: .seealso: NEPSetJacobian()
217: S*/
218: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(NEPJacobianFn)(NEP nep,PetscScalar lambda,Mat J,void *ctx);
220: SLEPC_EXTERN PetscErrorCode NEPSetFunction(NEP,Mat,Mat,NEPFunctionFn*,void*);
221: SLEPC_EXTERN PetscErrorCode NEPGetFunction(NEP,Mat*,Mat*,NEPFunctionFn**,void**);
222: SLEPC_EXTERN PetscErrorCode NEPSetJacobian(NEP,Mat,NEPJacobianFn*,void*);
223: SLEPC_EXTERN PetscErrorCode NEPGetJacobian(NEP,Mat*,NEPJacobianFn**,void**);
224: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "NEPSetFunction() and NEPSetJacobian()", ) static inline PetscErrorCode NEPSetDerivatives(NEP nep,Mat A,PetscErrorCode (*fun)(NEP,PetscScalar,PetscInt,Mat,void*),void *ctx) {(void)A;(void)fun;(void)ctx;SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented in this version");}
225: PETSC_DEPRECATED_FUNCTION(3, 12, 0, "NEPGetFunction() and NEPGetJacobian()", ) static inline PetscErrorCode NEPGetDerivatives(NEP nep,Mat *A,PetscErrorCode (**fun)(NEP,PetscScalar,PetscInt,Mat,void*),void **ctx) {(void)A;(void)fun;(void)ctx;SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Not implemented in this version");}
226: SLEPC_EXTERN PetscErrorCode NEPSetSplitOperator(NEP,PetscInt,Mat[],FN[],MatStructure);
227: SLEPC_EXTERN PetscErrorCode NEPGetSplitOperatorTerm(NEP,PetscInt,Mat*,FN*);
228: SLEPC_EXTERN PetscErrorCode NEPGetSplitOperatorInfo(NEP,PetscInt*,MatStructure*);
229: SLEPC_EXTERN PetscErrorCode NEPSetSplitPreconditioner(NEP,PetscInt,Mat[],MatStructure);
230: SLEPC_EXTERN PetscErrorCode NEPGetSplitPreconditionerTerm(NEP,PetscInt,Mat*);
231: SLEPC_EXTERN PetscErrorCode NEPGetSplitPreconditionerInfo(NEP,PetscInt*,MatStructure*);
233: SLEPC_EXTERN PetscErrorCode NEPSetBV(NEP,BV);
234: SLEPC_EXTERN PetscErrorCode NEPGetBV(NEP,BV*);
235: SLEPC_EXTERN PetscErrorCode NEPSetRG(NEP,RG);
236: SLEPC_EXTERN PetscErrorCode NEPGetRG(NEP,RG*);
237: SLEPC_EXTERN PetscErrorCode NEPSetDS(NEP,DS);
238: SLEPC_EXTERN PetscErrorCode NEPGetDS(NEP,DS*);
239: SLEPC_EXTERN PetscErrorCode NEPRefineGetKSP(NEP,KSP*);
240: SLEPC_EXTERN PetscErrorCode NEPSetTolerances(NEP,PetscReal,PetscInt);
241: SLEPC_EXTERN PetscErrorCode NEPGetTolerances(NEP,PetscReal*,PetscInt*);
242: SLEPC_EXTERN PetscErrorCode NEPSetConvergenceTest(NEP,NEPConv);
243: SLEPC_EXTERN PetscErrorCode NEPGetConvergenceTest(NEP,NEPConv*);
244: SLEPC_EXTERN PetscErrorCode NEPConvergedAbsolute(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
245: SLEPC_EXTERN PetscErrorCode NEPConvergedRelative(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
246: SLEPC_EXTERN PetscErrorCode NEPConvergedNorm(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
247: SLEPC_EXTERN PetscErrorCode NEPSetStoppingTest(NEP,NEPStop);
248: SLEPC_EXTERN PetscErrorCode NEPGetStoppingTest(NEP,NEPStop*);
249: SLEPC_EXTERN PetscErrorCode NEPStoppingBasic(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);
250: SLEPC_EXTERN PetscErrorCode NEPSetDimensions(NEP,PetscInt,PetscInt,PetscInt);
251: SLEPC_EXTERN PetscErrorCode NEPGetDimensions(NEP,PetscInt*,PetscInt*,PetscInt*);
252: SLEPC_EXTERN PetscErrorCode NEPSetRefine(NEP,NEPRefine,PetscInt,PetscReal,PetscInt,NEPRefineScheme);
253: SLEPC_EXTERN PetscErrorCode NEPGetRefine(NEP,NEPRefine*,PetscInt*,PetscReal*,PetscInt*,NEPRefineScheme*);
255: SLEPC_EXTERN PetscErrorCode NEPGetConverged(NEP,PetscInt*);
256: SLEPC_EXTERN PetscErrorCode NEPGetEigenpair(NEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
257: SLEPC_EXTERN PetscErrorCode NEPGetLeftEigenvector(NEP,PetscInt,Vec,Vec);
259: SLEPC_EXTERN PetscErrorCode NEPComputeError(NEP,PetscInt,NEPErrorType,PetscReal*);
260: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "NEPComputeError()", ) static inline PetscErrorCode NEPComputeRelativeError(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_RELATIVE,r);}
261: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "NEPComputeError() with NEP_ERROR_ABSOLUTE", ) static inline PetscErrorCode NEPComputeResidualNorm(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_ABSOLUTE,r);}
262: SLEPC_EXTERN PetscErrorCode NEPGetErrorEstimate(NEP,PetscInt,PetscReal*);
264: SLEPC_EXTERN PetscErrorCode NEPComputeFunction(NEP,PetscScalar,Mat,Mat);
265: SLEPC_EXTERN PetscErrorCode NEPComputeJacobian(NEP,PetscScalar,Mat);
266: SLEPC_EXTERN PetscErrorCode NEPApplyFunction(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
267: SLEPC_EXTERN PetscErrorCode NEPApplyAdjoint(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
268: SLEPC_EXTERN PetscErrorCode NEPApplyJacobian(NEP,PetscScalar,Vec,Vec,Vec,Mat);
269: SLEPC_EXTERN PetscErrorCode NEPProjectOperator(NEP,PetscInt,PetscInt);
270: SLEPC_EXTERN PetscErrorCode NEPGetIterationNumber(NEP,PetscInt*);
272: SLEPC_EXTERN PetscErrorCode NEPSetInitialSpace(NEP,PetscInt,Vec[]);
273: SLEPC_EXTERN PetscErrorCode NEPSetWhichEigenpairs(NEP,NEPWhich);
274: SLEPC_EXTERN PetscErrorCode NEPGetWhichEigenpairs(NEP,NEPWhich*);
275: SLEPC_EXTERN PetscErrorCode NEPSetTwoSided(NEP,PetscBool);
276: SLEPC_EXTERN PetscErrorCode NEPGetTwoSided(NEP,PetscBool*);
278: SLEPC_EXTERN PetscErrorCode NEPApplyResolvent(NEP,RG,PetscScalar,Vec,Vec);
280: SLEPC_EXTERN PetscErrorCode NEPSetTrackAll(NEP,PetscBool);
281: SLEPC_EXTERN PetscErrorCode NEPGetTrackAll(NEP,PetscBool*);
283: SLEPC_EXTERN PetscErrorCode NEPGetConvergedReason(NEP,NEPConvergedReason*);
285: SLEPC_EXTERN PetscErrorCode NEPMonitor(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
286: SLEPC_EXTERN PetscErrorCode NEPMonitorSet(NEP,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscCtxDestroyFn*);
287: SLEPC_EXTERN PetscErrorCode NEPMonitorCancel(NEP);
288: SLEPC_EXTERN PetscErrorCode NEPGetMonitorContext(NEP,void*);
290: SLEPC_EXTERN PetscErrorCode NEPMonitorSetFromOptions(NEP,const char[],const char[],void*,PetscBool);
291: SLEPC_EXTERN PetscErrorCode NEPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char*[],int,int,int,int,PetscDrawLG*);
292: SLEPC_EXTERN PetscErrorCode NEPMonitorFirst(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
293: SLEPC_EXTERN PetscErrorCode NEPMonitorFirstDrawLG(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
294: SLEPC_EXTERN PetscErrorCode NEPMonitorFirstDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
295: SLEPC_EXTERN PetscErrorCode NEPMonitorAll(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
296: SLEPC_EXTERN PetscErrorCode NEPMonitorAllDrawLG(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
297: SLEPC_EXTERN PetscErrorCode NEPMonitorAllDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
298: SLEPC_EXTERN PetscErrorCode NEPMonitorConverged(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
299: SLEPC_EXTERN PetscErrorCode NEPMonitorConvergedCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
300: SLEPC_EXTERN PetscErrorCode NEPMonitorConvergedDrawLG(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
301: SLEPC_EXTERN PetscErrorCode NEPMonitorConvergedDrawLGCreate(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
302: SLEPC_EXTERN PetscErrorCode NEPMonitorConvergedDestroy(PetscViewerAndFormat**);
304: SLEPC_EXTERN PetscErrorCode NEPSetOptionsPrefix(NEP,const char*);
305: SLEPC_EXTERN PetscErrorCode NEPAppendOptionsPrefix(NEP,const char*);
306: SLEPC_EXTERN PetscErrorCode NEPGetOptionsPrefix(NEP,const char*[]);
308: SLEPC_EXTERN PetscFunctionList NEPList;
309: SLEPC_EXTERN PetscFunctionList NEPMonitorList;
310: SLEPC_EXTERN PetscFunctionList NEPMonitorCreateList;
311: SLEPC_EXTERN PetscFunctionList NEPMonitorDestroyList;
312: SLEPC_EXTERN PetscErrorCode NEPRegister(const char[],PetscErrorCode(*)(NEP));
313: SLEPC_EXTERN PetscErrorCode NEPMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));
315: SLEPC_EXTERN PetscErrorCode NEPSetWorkVecs(NEP,PetscInt);
316: SLEPC_EXTERN PetscErrorCode NEPAllocateSolution(NEP,PetscInt);
318: /*S
319: NEPConvergenceTestFn - A prototype of a NEP convergence test function that would be passed to NEPSetConvergenceTestFunction()
321: Calling Sequence:
322: + nep - eigensolver context obtained from NEPCreate()
323: . eigr - real part of the eigenvalue
324: . eigi - imaginary part of the eigenvalue
325: . res - residual norm associated to the eigenpair
326: . errest - [output] computed error estimate
327: - ctx - [optional] user-defined context for private data for the
328: convergence test routine (may be NULL)
330: Level: advanced
332: .seealso: NEPSetConvergenceTestFunction()
333: S*/
334: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(NEPConvergenceTestFn)(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx);
336: /*S
337: NEPStoppingTestFn - A prototype of a NEP stopping test function that would be passed to NEPSetStoppingTestFunction()
339: Calling Sequence:
340: + nep - eigensolver context obtained from NEPCreate()
341: . its - current number of iterations
342: . max_it - maximum number of iterations
343: . nconv - number of currently converged eigenpairs
344: . nev - number of requested eigenpairs
345: . reason - [output] result of the stopping test
346: - ctx - [optional] user-defined context for private data for the
347: stopping test routine (may be NULL)
349: Level: advanced
351: .seealso: NEPSetStoppingTestFunction()
352: S*/
353: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(NEPStoppingTestFn)(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx);
355: SLEPC_EXTERN PetscErrorCode NEPSetConvergenceTestFunction(NEP,NEPConvergenceTestFn*,void*,PetscCtxDestroyFn*);
356: SLEPC_EXTERN PetscErrorCode NEPSetStoppingTestFunction(NEP,NEPStoppingTestFn*,void*,PetscCtxDestroyFn*);
357: SLEPC_EXTERN PetscErrorCode NEPSetEigenvalueComparison(NEP,SlepcEigenvalueComparisonFn*,void*);
359: /* --------- options specific to particular eigensolvers -------- */
361: SLEPC_EXTERN PetscErrorCode NEPRIISetMaximumIterations(NEP,PetscInt);
362: SLEPC_EXTERN PetscErrorCode NEPRIIGetMaximumIterations(NEP,PetscInt*);
363: SLEPC_EXTERN PetscErrorCode NEPRIISetLagPreconditioner(NEP,PetscInt);
364: SLEPC_EXTERN PetscErrorCode NEPRIIGetLagPreconditioner(NEP,PetscInt*);
365: SLEPC_EXTERN PetscErrorCode NEPRIISetConstCorrectionTol(NEP,PetscBool);
366: SLEPC_EXTERN PetscErrorCode NEPRIIGetConstCorrectionTol(NEP,PetscBool*);
367: SLEPC_EXTERN PetscErrorCode NEPRIISetHermitian(NEP,PetscBool);
368: SLEPC_EXTERN PetscErrorCode NEPRIIGetHermitian(NEP,PetscBool*);
369: SLEPC_EXTERN PetscErrorCode NEPRIISetDeflationThreshold(NEP,PetscReal);
370: SLEPC_EXTERN PetscErrorCode NEPRIIGetDeflationThreshold(NEP,PetscReal*);
371: SLEPC_EXTERN PetscErrorCode NEPRIISetKSP(NEP,KSP);
372: SLEPC_EXTERN PetscErrorCode NEPRIIGetKSP(NEP,KSP*);
374: SLEPC_EXTERN PetscErrorCode NEPSLPSetDeflationThreshold(NEP,PetscReal);
375: SLEPC_EXTERN PetscErrorCode NEPSLPGetDeflationThreshold(NEP,PetscReal*);
376: SLEPC_EXTERN PetscErrorCode NEPSLPSetEPS(NEP,EPS);
377: SLEPC_EXTERN PetscErrorCode NEPSLPGetEPS(NEP,EPS*);
378: SLEPC_EXTERN PetscErrorCode NEPSLPSetEPSLeft(NEP,EPS);
379: SLEPC_EXTERN PetscErrorCode NEPSLPGetEPSLeft(NEP,EPS*);
380: SLEPC_EXTERN PetscErrorCode NEPSLPSetKSP(NEP,KSP);
381: SLEPC_EXTERN PetscErrorCode NEPSLPGetKSP(NEP,KSP*);
383: SLEPC_EXTERN PetscErrorCode NEPNArnoldiSetKSP(NEP,KSP);
384: SLEPC_EXTERN PetscErrorCode NEPNArnoldiGetKSP(NEP,KSP*);
385: SLEPC_EXTERN PetscErrorCode NEPNArnoldiSetLagPreconditioner(NEP,PetscInt);
386: SLEPC_EXTERN PetscErrorCode NEPNArnoldiGetLagPreconditioner(NEP,PetscInt*);
388: /*E
389: NEPCISSExtraction - determines the extraction technique in the CISS solver
391: Level: advanced
393: .seealso: NEPCISSSetExtraction(), NEPCISSGetExtraction()
394: E*/
395: typedef enum { NEP_CISS_EXTRACTION_RITZ,
396: NEP_CISS_EXTRACTION_HANKEL,
397: NEP_CISS_EXTRACTION_CAA } NEPCISSExtraction;
398: SLEPC_EXTERN const char *NEPCISSExtractions[];
400: #if defined(PETSC_USE_COMPLEX) || defined(PETSC_CLANG_STATIC_ANALYZER)
401: SLEPC_EXTERN PetscErrorCode NEPCISSSetExtraction(NEP,NEPCISSExtraction);
402: SLEPC_EXTERN PetscErrorCode NEPCISSGetExtraction(NEP,NEPCISSExtraction*);
403: SLEPC_EXTERN PetscErrorCode NEPCISSSetSizes(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
404: SLEPC_EXTERN PetscErrorCode NEPCISSGetSizes(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
405: SLEPC_EXTERN PetscErrorCode NEPCISSSetThreshold(NEP,PetscReal,PetscReal);
406: SLEPC_EXTERN PetscErrorCode NEPCISSGetThreshold(NEP,PetscReal*,PetscReal*);
407: SLEPC_EXTERN PetscErrorCode NEPCISSSetRefinement(NEP,PetscInt,PetscInt);
408: SLEPC_EXTERN PetscErrorCode NEPCISSGetRefinement(NEP,PetscInt*,PetscInt*);
409: SLEPC_EXTERN PetscErrorCode NEPCISSGetKSPs(NEP,PetscInt*,KSP**);
410: #else
411: #define SlepcNEPCISSUnavailable(nep) do { \
412: PetscFunctionBegin; \
413: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"%s() not available with real scalars",PETSC_FUNCTION_NAME); \
414: } while (0)
415: static inline PetscErrorCode NEPCISSSetExtraction(NEP nep,PETSC_UNUSED NEPCISSExtraction ex) {SlepcNEPCISSUnavailable(nep);}
416: static inline PetscErrorCode NEPCISSGetExtraction(NEP nep,PETSC_UNUSED NEPCISSExtraction *ex) {SlepcNEPCISSUnavailable(nep);}
417: 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);}
418: 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);}
419: static inline PetscErrorCode NEPCISSSetThreshold(NEP nep,PETSC_UNUSED PetscReal delta,PETSC_UNUSED PetscReal spur) {SlepcNEPCISSUnavailable(nep);}
420: static inline PetscErrorCode NEPCISSGetThreshold(NEP nep,PETSC_UNUSED PetscReal *delta,PETSC_UNUSED PetscReal *spur) {SlepcNEPCISSUnavailable(nep);}
421: static inline PetscErrorCode NEPCISSSetRefinement(NEP nep,PETSC_UNUSED PetscInt inner,PETSC_UNUSED PetscInt blsize) {SlepcNEPCISSUnavailable(nep);}
422: static inline PetscErrorCode NEPCISSGetRefinement(NEP nep,PETSC_UNUSED PetscInt *inner,PETSC_UNUSED PetscInt *blsize) {SlepcNEPCISSUnavailable(nep);}
423: static inline PetscErrorCode NEPCISSGetKSPs(NEP nep,PETSC_UNUSED PetscInt *nsolve,PETSC_UNUSED KSP **ksp) {SlepcNEPCISSUnavailable(nep);}
424: #undef SlepcNEPCISSUnavailable
425: #endif
427: SLEPC_EXTERN PetscErrorCode NEPInterpolSetPEP(NEP,PEP);
428: SLEPC_EXTERN PetscErrorCode NEPInterpolGetPEP(NEP,PEP*);
429: SLEPC_EXTERN PetscErrorCode NEPInterpolSetInterpolation(NEP,PetscReal,PetscInt);
430: SLEPC_EXTERN PetscErrorCode NEPInterpolGetInterpolation(NEP,PetscReal*,PetscInt*);
432: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetSingularitiesFunction(NEP,PetscErrorCode (*)(NEP,PetscInt*,PetscScalar*,void*),void*);
433: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetSingularitiesFunction(NEP,PetscErrorCode (**)(NEP,PetscInt*,PetscScalar*,void*),void**);
434: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetRestart(NEP,PetscReal);
435: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetRestart(NEP,PetscReal*);
436: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetLocking(NEP,PetscBool);
437: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetLocking(NEP,PetscBool*);
438: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetInterpolation(NEP,PetscReal,PetscInt);
439: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetInterpolation(NEP,PetscReal*,PetscInt*);
440: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetRKShifts(NEP,PetscInt,PetscScalar[]);
441: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetRKShifts(NEP,PetscInt*,PetscScalar*[]);
442: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetKSPs(NEP,PetscInt*,KSP**);
443: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetFullBasis(NEP,PetscBool);
444: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetFullBasis(NEP,PetscBool*);
445: SLEPC_EXTERN PetscErrorCode NEPNLEIGSSetEPS(NEP,EPS);
446: SLEPC_EXTERN PetscErrorCode NEPNLEIGSGetEPS(NEP,EPS*);