Actual source code: slepcnep.h
slepc-3.6.1 2015-09-03
1: /*
2: User interface for SLEPc's nonlinear eigenvalue solvers.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
26: #include <slepceps.h>
27: #include <slepcpep.h>
28: #include <slepcfn.h>
30: PETSC_EXTERN PetscErrorCode NEPInitializePackage(void);
32: /*S
33: NEP - Abstract SLEPc object that manages all solvers for
34: nonlinear eigenvalue problems.
36: Level: beginner
38: .seealso: NEPCreate()
39: S*/
40: typedef struct _p_NEP* NEP;
42: /*J
43: NEPType - String with the name of a nonlinear eigensolver
45: Level: beginner
47: .seealso: NEPSetType(), NEP
48: J*/
49: typedef const char* NEPType;
50: #define NEPRII "rii"
51: #define NEPSLP "slp"
52: #define NEPNARNOLDI "narnoldi"
53: #define NEPCISS "ciss"
54: #define NEPINTERPOL "interpol"
56: /* Logging support */
57: PETSC_EXTERN PetscClassId NEP_CLASSID;
59: /*E
60: NEPWhich - Determines which part of the spectrum is requested
62: Level: intermediate
64: .seealso: NEPSetWhichEigenpairs(), NEPGetWhichEigenpairs()
65: E*/
66: typedef enum { NEP_LARGEST_MAGNITUDE=1,
67: NEP_SMALLEST_MAGNITUDE,
68: NEP_LARGEST_REAL,
69: NEP_SMALLEST_REAL,
70: NEP_LARGEST_IMAGINARY,
71: NEP_SMALLEST_IMAGINARY,
72: NEP_TARGET_MAGNITUDE,
73: NEP_TARGET_REAL,
74: NEP_TARGET_IMAGINARY} NEPWhich;
76: /*E
77: NEPErrorType - The error type used to assess accuracy of computed solutions
79: Level: intermediate
81: .seealso: NEPComputeError()
82: E*/
83: typedef enum { NEP_ERROR_ABSOLUTE,
84: NEP_ERROR_RELATIVE } NEPErrorType;
85: PETSC_EXTERN const char *NEPErrorTypes[];
87: /*E
88: NEPRefine - The refinement type
90: Level: intermediate
92: .seealso: NEPSetRefine()
93: E*/
94: typedef enum { NEP_REFINE_NONE,
95: NEP_REFINE_SIMPLE,
96: NEP_REFINE_MULTIPLE } NEPRefine;
97: PETSC_EXTERN const char *NEPRefineTypes[];
99: /*E
100: NEPConvergedReason - Reason a nonlinear eigensolver was said to
101: have converged or diverged
103: Level: intermediate
105: .seealso: NEPSolve(), NEPGetConvergedReason(), NEPSetTolerances()
106: E*/
107: typedef enum {/* converged */
108: NEP_CONVERGED_FNORM_ABS = 2,
109: NEP_CONVERGED_FNORM_RELATIVE = 3,
110: NEP_CONVERGED_SNORM_RELATIVE = 4,
111: /* diverged */
112: NEP_DIVERGED_LINEAR_SOLVE = -1,
113: NEP_DIVERGED_FUNCTION_COUNT = -2,
114: NEP_DIVERGED_MAX_IT = -3,
115: NEP_DIVERGED_BREAKDOWN = -4,
116: NEP_DIVERGED_FNORM_NAN = -5,
117: NEP_CONVERGED_ITERATING = 0} NEPConvergedReason;
118: PETSC_EXTERN const char *const*NEPConvergedReasons;
120: PETSC_EXTERN PetscErrorCode NEPCreate(MPI_Comm,NEP*);
121: PETSC_EXTERN PetscErrorCode NEPDestroy(NEP*);
122: PETSC_EXTERN PetscErrorCode NEPReset(NEP);
123: PETSC_EXTERN PetscErrorCode NEPSetType(NEP,NEPType);
124: PETSC_EXTERN PetscErrorCode NEPGetType(NEP,NEPType*);
125: PETSC_EXTERN PetscErrorCode NEPSetTarget(NEP,PetscScalar);
126: PETSC_EXTERN PetscErrorCode NEPGetTarget(NEP,PetscScalar*);
127: PETSC_EXTERN PetscErrorCode NEPSetKSP(NEP,KSP);
128: PETSC_EXTERN PetscErrorCode NEPGetKSP(NEP,KSP*);
129: PETSC_EXTERN PetscErrorCode NEPSetFromOptions(NEP);
130: PETSC_EXTERN PetscErrorCode NEPSetUp(NEP);
131: PETSC_EXTERN PetscErrorCode NEPSolve(NEP);
132: PETSC_EXTERN PetscErrorCode NEPView(NEP,PetscViewer);
133: PETSC_STATIC_INLINE PetscErrorCode NEPViewFromOptions(NEP nep,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)nep,obj,name);}
134: PETSC_EXTERN PetscErrorCode NEPErrorView(NEP,NEPErrorType,PetscViewer);
135: PETSC_EXTERN PetscErrorCode NEPErrorViewFromOptions(NEP);
136: PETSC_EXTERN PetscErrorCode NEPReasonView(NEP,PetscViewer);
137: PETSC_EXTERN PetscErrorCode NEPReasonViewFromOptions(NEP);
138: PETSC_EXTERN PetscErrorCode NEPValuesView(NEP,PetscViewer);
139: PETSC_EXTERN PetscErrorCode NEPValuesViewFromOptions(NEP);
140: PETSC_EXTERN PetscErrorCode NEPVectorsView(NEP,PetscViewer);
141: PETSC_EXTERN PetscErrorCode NEPVectorsViewFromOptions(NEP);
143: PETSC_EXTERN PetscErrorCode NEPSetFunction(NEP,Mat,Mat,PetscErrorCode (*)(NEP,PetscScalar,Mat,Mat,void*),void*);
144: PETSC_EXTERN PetscErrorCode NEPGetFunction(NEP,Mat*,Mat*,PetscErrorCode (**)(NEP,PetscScalar,Mat,Mat,void*),void**);
145: PETSC_EXTERN PetscErrorCode NEPSetJacobian(NEP,Mat,PetscErrorCode (*)(NEP,PetscScalar,Mat,void*),void*);
146: PETSC_EXTERN PetscErrorCode NEPGetJacobian(NEP,Mat*,PetscErrorCode (**)(NEP,PetscScalar,Mat,void*),void**);
147: PETSC_EXTERN PetscErrorCode NEPSetSplitOperator(NEP,PetscInt,Mat*,FN*,MatStructure);
148: PETSC_EXTERN PetscErrorCode NEPGetSplitOperatorTerm(NEP,PetscInt,Mat*,FN*);
149: PETSC_EXTERN PetscErrorCode NEPGetSplitOperatorInfo(NEP,PetscInt*,MatStructure*);
151: PETSC_EXTERN PetscErrorCode NEPSetBV(NEP,BV);
152: PETSC_EXTERN PetscErrorCode NEPGetBV(NEP,BV*);
153: PETSC_EXTERN PetscErrorCode NEPSetRG(NEP,RG);
154: PETSC_EXTERN PetscErrorCode NEPGetRG(NEP,RG*);
155: PETSC_EXTERN PetscErrorCode NEPSetDS(NEP,DS);
156: PETSC_EXTERN PetscErrorCode NEPGetDS(NEP,DS*);
157: PETSC_EXTERN PetscErrorCode NEPSetTolerances(NEP,PetscReal,PetscReal,PetscReal,PetscInt,PetscInt);
158: PETSC_EXTERN PetscErrorCode NEPGetTolerances(NEP,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt*);
159: PETSC_EXTERN PetscErrorCode NEPSetConvergenceTest(NEP,PetscErrorCode (*)(NEP,PetscInt,PetscReal,PetscReal,PetscReal,NEPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
160: PETSC_EXTERN PetscErrorCode NEPConvergedDefault(NEP,PetscInt,PetscReal,PetscReal,PetscReal,NEPConvergedReason*,void*);
161: PETSC_EXTERN PetscErrorCode NEPSetDimensions(NEP,PetscInt,PetscInt,PetscInt);
162: PETSC_EXTERN PetscErrorCode NEPGetDimensions(NEP,PetscInt*,PetscInt*,PetscInt*);
163: PETSC_EXTERN PetscErrorCode NEPSetRefine(NEP,NEPRefine,PetscInt,PetscReal,PetscInt);
164: PETSC_EXTERN PetscErrorCode NEPGetRefine(NEP,NEPRefine*,PetscInt*,PetscReal*,PetscInt*);
165: PETSC_EXTERN PetscErrorCode NEPSetLagPreconditioner(NEP,PetscInt);
166: PETSC_EXTERN PetscErrorCode NEPGetLagPreconditioner(NEP,PetscInt*);
167: PETSC_EXTERN PetscErrorCode NEPSetConstCorrectionTol(NEP,PetscBool);
168: PETSC_EXTERN PetscErrorCode NEPGetConstCorrectionTol(NEP,PetscBool*);
170: PETSC_EXTERN PetscErrorCode NEPGetConverged(NEP,PetscInt*);
171: PETSC_EXTERN PetscErrorCode NEPGetEigenpair(NEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
173: PETSC_EXTERN PetscErrorCode NEPComputeError(NEP,PetscInt,NEPErrorType,PetscReal*);
174: PETSC_DEPRECATED("Use NEPComputeError()") PETSC_STATIC_INLINE PetscErrorCode NEPComputeRelativeError(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_RELATIVE,r);}
175: PETSC_DEPRECATED("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);}
176: PETSC_EXTERN PetscErrorCode NEPGetErrorEstimate(NEP,PetscInt,PetscReal*);
178: PETSC_EXTERN PetscErrorCode NEPComputeFunction(NEP,PetscScalar,Mat,Mat);
179: PETSC_EXTERN PetscErrorCode NEPComputeJacobian(NEP,PetscScalar,Mat);
180: PETSC_EXTERN PetscErrorCode NEPApplyFunction(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
181: PETSC_EXTERN PetscErrorCode NEPApplyJacobian(NEP,PetscScalar,Vec,Vec,Vec,Mat);
182: PETSC_EXTERN PetscErrorCode NEPProjectOperator(NEP,PetscInt,PetscInt);
184: PETSC_EXTERN PetscErrorCode NEPMonitor(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt);
185: PETSC_EXTERN PetscErrorCode NEPMonitorSet(NEP,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
186: PETSC_EXTERN PetscErrorCode NEPMonitorCancel(NEP);
187: PETSC_EXTERN PetscErrorCode NEPGetMonitorContext(NEP,void **);
188: PETSC_EXTERN PetscErrorCode NEPGetIterationNumber(NEP,PetscInt*);
190: PETSC_EXTERN PetscErrorCode NEPSetInitialSpace(NEP,PetscInt,Vec*);
191: PETSC_EXTERN PetscErrorCode NEPSetWhichEigenpairs(NEP,NEPWhich);
192: PETSC_EXTERN PetscErrorCode NEPGetWhichEigenpairs(NEP,NEPWhich*);
194: PETSC_EXTERN PetscErrorCode NEPMonitorAll(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt,void*);
195: PETSC_EXTERN PetscErrorCode NEPMonitorFirst(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt,void*);
196: PETSC_EXTERN PetscErrorCode NEPMonitorConverged(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt,void*);
197: PETSC_EXTERN PetscErrorCode NEPMonitorLG(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt,void*);
198: PETSC_EXTERN PetscErrorCode NEPMonitorLGAll(NEP,PetscInt,PetscInt,PetscScalar*,PetscReal*,PetscInt,void*);
200: PETSC_EXTERN PetscErrorCode NEPSetTrackAll(NEP,PetscBool);
201: PETSC_EXTERN PetscErrorCode NEPGetTrackAll(NEP,PetscBool*);
203: PETSC_EXTERN PetscErrorCode NEPSetOptionsPrefix(NEP,const char*);
204: PETSC_EXTERN PetscErrorCode NEPAppendOptionsPrefix(NEP,const char*);
205: PETSC_EXTERN PetscErrorCode NEPGetOptionsPrefix(NEP,const char*[]);
207: PETSC_EXTERN PetscErrorCode NEPGetConvergedReason(NEP,NEPConvergedReason *);
209: PETSC_EXTERN PetscFunctionList NEPList;
210: PETSC_EXTERN PetscErrorCode NEPRegister(const char[],PetscErrorCode(*)(NEP));
212: PETSC_EXTERN PetscErrorCode NEPSetWorkVecs(NEP,PetscInt);
213: PETSC_EXTERN PetscErrorCode NEPAllocateSolution(NEP,PetscInt);
215: /* --------- options specific to particular eigensolvers -------- */
217: PETSC_EXTERN PetscErrorCode NEPSLPSetEPS(NEP,EPS);
218: PETSC_EXTERN PetscErrorCode NEPSLPGetEPS(NEP,EPS*);
220: PETSC_EXTERN PetscErrorCode NEPCISSSetSizes(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
221: PETSC_EXTERN PetscErrorCode NEPCISSGetSizes(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
222: PETSC_EXTERN PetscErrorCode NEPCISSSetThreshold(NEP,PetscReal,PetscReal);
223: PETSC_EXTERN PetscErrorCode NEPCISSGetThreshold(NEP,PetscReal*,PetscReal*);
224: PETSC_EXTERN PetscErrorCode NEPCISSSetRefinement(NEP,PetscInt,PetscInt,PetscInt);
225: PETSC_EXTERN PetscErrorCode NEPCISSGetRefinement(NEP,PetscInt*,PetscInt*,PetscInt*);
227: PETSC_EXTERN PetscErrorCode NEPInterpolSetPEP(NEP,PEP);
228: PETSC_EXTERN PetscErrorCode NEPInterpolGetPEP(NEP,PEP*);
229: PETSC_EXTERN PetscErrorCode NEPInterpolSetDegree(NEP,PetscInt);
230: PETSC_EXTERN PetscErrorCode NEPInterpolGetDegree(NEP,PetscInt*);
232: #endif