Actual source code: slepcnep.h

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    User interface for SLEPc's nonlinear eigenvalue solvers.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, 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"
 55: #define NEPNLEIGS    "nleigs"

 57: /* Logging support */
 58: PETSC_EXTERN PetscClassId NEP_CLASSID;

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

 63:     Level: intermediate

 65: .seealso: NEPSetWhichEigenpairs(), NEPGetWhichEigenpairs()
 66: E*/
 67: typedef enum { NEP_LARGEST_MAGNITUDE=1,
 68:                NEP_SMALLEST_MAGNITUDE,
 69:                NEP_LARGEST_REAL,
 70:                NEP_SMALLEST_REAL,
 71:                NEP_LARGEST_IMAGINARY,
 72:                NEP_SMALLEST_IMAGINARY,
 73:                NEP_TARGET_MAGNITUDE,
 74:                NEP_TARGET_REAL,
 75:                NEP_TARGET_IMAGINARY,
 76:                NEP_ALL,
 77:                NEP_WHICH_USER } NEPWhich;

 79: /*E
 80:     NEPErrorType - The error type used to assess accuracy of computed solutions

 82:     Level: intermediate

 84: .seealso: NEPComputeError()
 85: E*/
 86: typedef enum { NEP_ERROR_ABSOLUTE,
 87:                NEP_ERROR_RELATIVE,
 88:                NEP_ERROR_BACKWARD } NEPErrorType;
 89: PETSC_EXTERN const char *NEPErrorTypes[];

 91: /*E
 92:     NEPRefine - The refinement type

 94:     Level: intermediate

 96: .seealso: NEPSetRefine()
 97: E*/
 98: typedef enum { NEP_REFINE_NONE,
 99:                NEP_REFINE_SIMPLE,
100:                NEP_REFINE_MULTIPLE } NEPRefine;
101: PETSC_EXTERN const char *NEPRefineTypes[];

103: /*E
104:     NEPRefineScheme - The scheme used for solving linear systems during iterative refinement

106:     Level: intermediate

108: .seealso: NEPSetRefine()
109: E*/
110: typedef enum { NEP_REFINE_SCHEME_SCHUR=1,
111:                NEP_REFINE_SCHEME_MBE,
112:                NEP_REFINE_SCHEME_EXPLICIT } NEPRefineScheme;
113: PETSC_EXTERN const char *NEPRefineSchemes[];

115: /*E
116:     NEPConv - Determines the convergence test

118:     Level: intermediate

120: .seealso: NEPSetConvergenceTest(), NEPSetConvergenceTestFunction()
121: E*/
122: typedef enum { NEP_CONV_ABS,
123:                NEP_CONV_REL,
124:                NEP_CONV_NORM,
125:                NEP_CONV_USER } NEPConv;

127: /*E
128:     NEPStop - Determines the stopping test

130:     Level: advanced

132: .seealso: NEPSetStoppingTest(), NEPSetStoppingTestFunction()
133: E*/
134: typedef enum { NEP_STOP_BASIC,
135:                NEP_STOP_USER } NEPStop;

137: /*E
138:     NEPConvergedReason - Reason a nonlinear eigensolver was said to
139:          have converged or diverged

141:     Level: intermediate

143: .seealso: NEPSolve(), NEPGetConvergedReason(), NEPSetTolerances()
144: E*/
145: typedef enum {/* converged */
146:               NEP_CONVERGED_TOL                =  1,
147:               NEP_CONVERGED_USER               =  2,
148:               /* diverged */
149:               NEP_DIVERGED_ITS                 = -1,
150:               NEP_DIVERGED_BREAKDOWN           = -2,
151:                     /* unused                  = -3 */
152:               NEP_DIVERGED_LINEAR_SOLVE        = -4,
153:               NEP_CONVERGED_ITERATING          =  0} NEPConvergedReason;
154: PETSC_EXTERN const char *const*NEPConvergedReasons;

156: PETSC_EXTERN PetscErrorCode NEPCreate(MPI_Comm,NEP*);
157: PETSC_EXTERN PetscErrorCode NEPDestroy(NEP*);
158: PETSC_EXTERN PetscErrorCode NEPReset(NEP);
159: PETSC_EXTERN PetscErrorCode NEPSetType(NEP,NEPType);
160: PETSC_EXTERN PetscErrorCode NEPGetType(NEP,NEPType*);
161: PETSC_EXTERN PetscErrorCode NEPSetTarget(NEP,PetscScalar);
162: PETSC_EXTERN PetscErrorCode NEPGetTarget(NEP,PetscScalar*);
163: PETSC_EXTERN PetscErrorCode NEPSetFromOptions(NEP);
164: PETSC_EXTERN PetscErrorCode NEPSetUp(NEP);
165: PETSC_EXTERN PetscErrorCode NEPSolve(NEP);
166: PETSC_EXTERN PetscErrorCode NEPView(NEP,PetscViewer);
167: PETSC_STATIC_INLINE PetscErrorCode NEPViewFromOptions(NEP nep,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)nep,obj,name);}
168: PETSC_EXTERN PetscErrorCode NEPErrorView(NEP,NEPErrorType,PetscViewer);
169: PETSC_EXTERN PetscErrorCode NEPErrorViewFromOptions(NEP);
170: PETSC_EXTERN PetscErrorCode NEPReasonView(NEP,PetscViewer);
171: PETSC_EXTERN PetscErrorCode NEPReasonViewFromOptions(NEP);
172: PETSC_EXTERN PetscErrorCode NEPValuesView(NEP,PetscViewer);
173: PETSC_EXTERN PetscErrorCode NEPValuesViewFromOptions(NEP);
174: PETSC_EXTERN PetscErrorCode NEPVectorsView(NEP,PetscViewer);
175: PETSC_EXTERN PetscErrorCode NEPVectorsViewFromOptions(NEP);

177: PETSC_EXTERN PetscErrorCode NEPSetFunction(NEP,Mat,Mat,PetscErrorCode (*)(NEP,PetscScalar,Mat,Mat,void*),void*);
178: PETSC_EXTERN PetscErrorCode NEPGetFunction(NEP,Mat*,Mat*,PetscErrorCode (**)(NEP,PetscScalar,Mat,Mat,void*),void**);
179: PETSC_EXTERN PetscErrorCode NEPSetJacobian(NEP,Mat,PetscErrorCode (*)(NEP,PetscScalar,Mat,void*),void*);
180: PETSC_EXTERN PetscErrorCode NEPGetJacobian(NEP,Mat*,PetscErrorCode (**)(NEP,PetscScalar,Mat,void*),void**);
181: PETSC_EXTERN PetscErrorCode NEPSetDerivatives(NEP,Mat,PetscErrorCode (*)(NEP,PetscScalar,PetscInt,Mat,void*),void*);
182: PETSC_EXTERN PetscErrorCode NEPGetDerivatives(NEP,Mat*,PetscErrorCode (**)(NEP,PetscScalar,PetscInt,Mat,void*),void**);
183: PETSC_EXTERN PetscErrorCode NEPSetSplitOperator(NEP,PetscInt,Mat*,FN*,MatStructure);
184: PETSC_EXTERN PetscErrorCode NEPGetSplitOperatorTerm(NEP,PetscInt,Mat*,FN*);
185: PETSC_EXTERN PetscErrorCode NEPGetSplitOperatorInfo(NEP,PetscInt*,MatStructure*);

187: PETSC_EXTERN PetscErrorCode NEPSetBV(NEP,BV);
188: PETSC_EXTERN PetscErrorCode NEPGetBV(NEP,BV*);
189: PETSC_EXTERN PetscErrorCode NEPSetRG(NEP,RG);
190: PETSC_EXTERN PetscErrorCode NEPGetRG(NEP,RG*);
191: PETSC_EXTERN PetscErrorCode NEPSetDS(NEP,DS);
192: PETSC_EXTERN PetscErrorCode NEPGetDS(NEP,DS*);
193: PETSC_EXTERN PetscErrorCode NEPRefineGetKSP(NEP,KSP*);
194: PETSC_EXTERN PetscErrorCode NEPSetTolerances(NEP,PetscReal,PetscInt);
195: PETSC_EXTERN PetscErrorCode NEPGetTolerances(NEP,PetscReal*,PetscInt*);
196: PETSC_EXTERN PetscErrorCode NEPSetConvergenceTestFunction(NEP,PetscErrorCode (*)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void*,PetscErrorCode (*)(void*));
197: PETSC_EXTERN PetscErrorCode NEPSetConvergenceTest(NEP,NEPConv);
198: PETSC_EXTERN PetscErrorCode NEPGetConvergenceTest(NEP,NEPConv*);
199: PETSC_EXTERN PetscErrorCode NEPConvergedAbsolute(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
200: PETSC_EXTERN PetscErrorCode NEPConvergedRelative(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
201: PETSC_EXTERN PetscErrorCode NEPConvergedNorm(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*);
202: PETSC_EXTERN PetscErrorCode NEPSetStoppingTestFunction(NEP,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
203: PETSC_EXTERN PetscErrorCode NEPSetStoppingTest(NEP,NEPStop);
204: PETSC_EXTERN PetscErrorCode NEPGetStoppingTest(NEP,NEPStop*);
205: PETSC_EXTERN PetscErrorCode NEPStoppingBasic(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*);
206: PETSC_EXTERN PetscErrorCode NEPSetDimensions(NEP,PetscInt,PetscInt,PetscInt);
207: PETSC_EXTERN PetscErrorCode NEPGetDimensions(NEP,PetscInt*,PetscInt*,PetscInt*);
208: PETSC_EXTERN PetscErrorCode NEPSetRefine(NEP,NEPRefine,PetscInt,PetscReal,PetscInt,NEPRefineScheme);
209: PETSC_EXTERN PetscErrorCode NEPGetRefine(NEP,NEPRefine*,PetscInt*,PetscReal*,PetscInt*,NEPRefineScheme*);

211: PETSC_EXTERN PetscErrorCode NEPGetConverged(NEP,PetscInt*);
212: PETSC_EXTERN PetscErrorCode NEPGetEigenpair(NEP,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);

214: PETSC_EXTERN PetscErrorCode NEPComputeError(NEP,PetscInt,NEPErrorType,PetscReal*);
215: PETSC_DEPRECATED("Use NEPComputeError()") PETSC_STATIC_INLINE PetscErrorCode NEPComputeRelativeError(NEP nep,PetscInt i,PetscReal *r) {return NEPComputeError(nep,i,NEP_ERROR_RELATIVE,r);}
216: 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);}
217: PETSC_EXTERN PetscErrorCode NEPGetErrorEstimate(NEP,PetscInt,PetscReal*);

219: PETSC_EXTERN PetscErrorCode NEPComputeFunction(NEP,PetscScalar,Mat,Mat);
220: PETSC_EXTERN PetscErrorCode NEPComputeJacobian(NEP,PetscScalar,Mat);
221: PETSC_EXTERN PetscErrorCode NEPApplyFunction(NEP,PetscScalar,Vec,Vec,Vec,Mat,Mat);
222: PETSC_EXTERN PetscErrorCode NEPApplyJacobian(NEP,PetscScalar,Vec,Vec,Vec,Mat);
223: PETSC_EXTERN PetscErrorCode NEPProjectOperator(NEP,PetscInt,PetscInt);

225: PETSC_EXTERN PetscErrorCode NEPMonitor(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt);
226: PETSC_EXTERN PetscErrorCode NEPMonitorSet(NEP,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void*,PetscErrorCode (*)(void**));
227: PETSC_EXTERN PetscErrorCode NEPMonitorSetFromOptions(NEP,const char*,const char*,const char*,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool);
228: PETSC_EXTERN PetscErrorCode NEPConvMonitorSetFromOptions(NEP,const char*,const char*,const char*,PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor));
229: PETSC_EXTERN PetscErrorCode NEPMonitorCancel(NEP);
230: PETSC_EXTERN PetscErrorCode NEPGetMonitorContext(NEP,void **);
231: PETSC_EXTERN PetscErrorCode NEPGetIterationNumber(NEP,PetscInt*);

233: PETSC_EXTERN PetscErrorCode NEPSetInitialSpace(NEP,PetscInt,Vec*);
234: PETSC_EXTERN PetscErrorCode NEPSetWhichEigenpairs(NEP,NEPWhich);
235: PETSC_EXTERN PetscErrorCode NEPGetWhichEigenpairs(NEP,NEPWhich*);
236: PETSC_EXTERN PetscErrorCode NEPSetEigenvalueComparison(NEP,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void*);

238: PETSC_EXTERN PetscErrorCode NEPMonitorAll(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
239: PETSC_EXTERN PetscErrorCode NEPMonitorFirst(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*);
240: PETSC_EXTERN PetscErrorCode NEPMonitorConverged(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor);
241: PETSC_EXTERN PetscErrorCode NEPMonitorLGCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*);
242: PETSC_EXTERN PetscErrorCode NEPMonitorLG(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
243: PETSC_EXTERN PetscErrorCode NEPMonitorLGAll(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);

245: PETSC_EXTERN PetscErrorCode NEPSetTrackAll(NEP,PetscBool);
246: PETSC_EXTERN PetscErrorCode NEPGetTrackAll(NEP,PetscBool*);

248: PETSC_EXTERN PetscErrorCode NEPSetOptionsPrefix(NEP,const char*);
249: PETSC_EXTERN PetscErrorCode NEPAppendOptionsPrefix(NEP,const char*);
250: PETSC_EXTERN PetscErrorCode NEPGetOptionsPrefix(NEP,const char*[]);

252: PETSC_EXTERN PetscErrorCode NEPGetConvergedReason(NEP,NEPConvergedReason *);

254: PETSC_EXTERN PetscFunctionList NEPList;
255: PETSC_EXTERN PetscErrorCode NEPRegister(const char[],PetscErrorCode(*)(NEP));

257: PETSC_EXTERN PetscErrorCode NEPSetWorkVecs(NEP,PetscInt);
258: PETSC_EXTERN PetscErrorCode NEPAllocateSolution(NEP,PetscInt);

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

262: PETSC_EXTERN PetscErrorCode NEPRIISetMaximumIterations(NEP,PetscInt);
263: PETSC_EXTERN PetscErrorCode NEPRIIGetMaximumIterations(NEP,PetscInt*);
264: PETSC_EXTERN PetscErrorCode NEPRIISetLagPreconditioner(NEP,PetscInt);
265: PETSC_EXTERN PetscErrorCode NEPRIIGetLagPreconditioner(NEP,PetscInt*);
266: PETSC_EXTERN PetscErrorCode NEPRIISetConstCorrectionTol(NEP,PetscBool);
267: PETSC_EXTERN PetscErrorCode NEPRIIGetConstCorrectionTol(NEP,PetscBool*);
268: PETSC_EXTERN PetscErrorCode NEPRIISetKSP(NEP,KSP);
269: PETSC_EXTERN PetscErrorCode NEPRIIGetKSP(NEP,KSP*);

271: PETSC_EXTERN PetscErrorCode NEPSLPSetEPS(NEP,EPS);
272: PETSC_EXTERN PetscErrorCode NEPSLPGetEPS(NEP,EPS*);

274: PETSC_EXTERN PetscErrorCode NEPNArnoldiSetKSP(NEP,KSP);
275: PETSC_EXTERN PetscErrorCode NEPNArnoldiGetKSP(NEP,KSP*);

277: PETSC_EXTERN PetscErrorCode NEPCISSSetSizes(NEP,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
278: PETSC_EXTERN PetscErrorCode NEPCISSGetSizes(NEP,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
279: PETSC_EXTERN PetscErrorCode NEPCISSSetThreshold(NEP,PetscReal,PetscReal);
280: PETSC_EXTERN PetscErrorCode NEPCISSGetThreshold(NEP,PetscReal*,PetscReal*);
281: PETSC_EXTERN PetscErrorCode NEPCISSSetRefinement(NEP,PetscInt,PetscInt);
282: PETSC_EXTERN PetscErrorCode NEPCISSGetRefinement(NEP,PetscInt*,PetscInt*);

284: PETSC_EXTERN PetscErrorCode NEPInterpolSetPEP(NEP,PEP);
285: PETSC_EXTERN PetscErrorCode NEPInterpolGetPEP(NEP,PEP*);
286: PETSC_EXTERN PetscErrorCode NEPInterpolSetDegree(NEP,PetscInt);
287: PETSC_EXTERN PetscErrorCode NEPInterpolGetDegree(NEP,PetscInt*);

289: PETSC_EXTERN PetscErrorCode NEPNLEIGSSetSingularitiesFunction(NEP,PetscErrorCode (*)(NEP,PetscInt*,PetscScalar*,void*),void*);
290: PETSC_EXTERN PetscErrorCode NEPNLEIGSGetSingularitiesFunction(NEP,PetscErrorCode (**)(NEP,PetscInt*,PetscScalar*,void*),void **);
291: PETSC_EXTERN PetscErrorCode NEPNLEIGSSetRestart(NEP,PetscReal);
292: PETSC_EXTERN PetscErrorCode NEPNLEIGSGetRestart(NEP,PetscReal*);
293: PETSC_EXTERN PetscErrorCode NEPNLEIGSSetLocking(NEP,PetscBool);
294: PETSC_EXTERN PetscErrorCode NEPNLEIGSGetLocking(NEP,PetscBool*);
295: PETSC_EXTERN PetscErrorCode NEPNLEIGSSetInterpolation(NEP,PetscReal,PetscInt);
296: PETSC_EXTERN PetscErrorCode NEPNLEIGSGetInterpolation(NEP,PetscReal*,PetscInt*);
297: PETSC_EXTERN PetscErrorCode NEPNLEIGSSetTrueResidual(NEP,PetscBool);
298: PETSC_EXTERN PetscErrorCode NEPNLEIGSGetTrueResidual(NEP,PetscBool*);
299: PETSC_EXTERN PetscErrorCode NEPNLEIGSSetRKShifts(NEP,PetscInt,PetscScalar*);
300: PETSC_EXTERN PetscErrorCode NEPNLEIGSGetRKShifts(NEP,PetscInt*,PetscScalar**);
301: PETSC_EXTERN PetscErrorCode NEPNLEIGSGetKSPs(NEP,KSP**);

303: #endif