Actual source code: petscsnes.h

  1: /*
  2:     User interface for the nonlinear solvers package.
  3: */
 6:  #include petscksp.h
  7: PETSC_EXTERN_CXX_BEGIN

  9: /*S
 10:      SNES - Abstract PETSc object that manages all nonlinear solves

 12:    Level: beginner

 14:   Concepts: nonlinear solvers

 16: .seealso:  SNESCreate(), SNESSetType(), SNESType, TS, KSP, KSP, PC
 17: S*/
 18: typedef struct _p_SNES* SNES;

 20: /*J
 21:     SNESType - String with the name of a PETSc SNES method or the creation function
 22:        with an optional dynamic library name, for example
 23:        http://www.mcs.anl.gov/petsc/lib.a:mysnescreate()

 25:    Level: beginner

 27: .seealso: SNESSetType(), SNES
 28: J*/
 29: #define SNESType char*
 30: #define SNESLS           "ls"
 31: #define SNESTR           "tr"
 32: #define SNESPYTHON       "python"
 33: #define SNESTEST         "test"
 34: #define SNESNRICHARDSON  "nrichardson"
 35: #define SNESKSPONLY      "ksponly"
 36: #define SNESVIRS         "virs"
 37: #define SNESVISS         "viss"
 38: #define SNESNGMRES       "ngmres"
 39: #define SNESQN           "qn"
 40: #define SNESSHELL        "shell"
 41: #define SNESGS           "gs"
 42: #define SNESNCG          "ncg"
 43: #define SNESSORQN        "sorqn"
 44: #define SNESFAS          "fas"
 45: #define SNESMS           "ms"

 47: /* Logging support */
 48: extern PetscClassId  SNES_CLASSID;

 50: extern PetscErrorCode  SNESInitializePackage(const char[]);

 52: extern PetscErrorCode  SNESCreate(MPI_Comm,SNES*);
 53: extern PetscErrorCode  SNESReset(SNES);
 54: extern PetscErrorCode  SNESDestroy(SNES*);
 55: extern PetscErrorCode  SNESSetType(SNES,const SNESType);
 56: extern PetscErrorCode  SNESMonitor(SNES,PetscInt,PetscReal);
 57: extern PetscErrorCode  SNESMonitorSet(SNES,PetscErrorCode(*)(SNES,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**));
 58: extern PetscErrorCode  SNESMonitorCancel(SNES);
 59: extern PetscErrorCode  SNESSetConvergenceHistory(SNES,PetscReal[],PetscInt[],PetscInt,PetscBool );
 60: extern PetscErrorCode  SNESGetConvergenceHistory(SNES,PetscReal*[],PetscInt *[],PetscInt *);
 61: extern PetscErrorCode  SNESSetUp(SNES);
 62: extern PetscErrorCode  SNESSolve(SNES,Vec,Vec);
 63: extern PetscErrorCode  SNESSetErrorIfNotConverged(SNES,PetscBool );
 64: extern PetscErrorCode  SNESGetErrorIfNotConverged(SNES,PetscBool  *);


 67: extern PetscErrorCode  SNESAddOptionsChecker(PetscErrorCode (*)(SNES));

 69: extern PetscErrorCode  SNESSetUpdate(SNES, PetscErrorCode (*)(SNES, PetscInt));
 70: extern PetscErrorCode  SNESDefaultUpdate(SNES, PetscInt);

 72: extern PetscFList SNESList;
 73: extern PetscErrorCode  SNESRegisterDestroy(void);
 74: extern PetscErrorCode  SNESRegisterAll(const char[]);

 76: extern PetscErrorCode  SNESRegister(const char[],const char[],const char[],PetscErrorCode (*)(SNES));

 78: /*MC
 79:    SNESRegisterDynamic - Adds a method to the nonlinear solver package.

 81:    Synopsis:
 82:    PetscErrorCode SNESRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(SNES))

 84:    Not collective

 86:    Input Parameters:
 87: +  name_solver - name of a new user-defined solver
 88: .  path - path (either absolute or relative) the library containing this solver
 89: .  name_create - name of routine to create method context
 90: -  routine_create - routine to create method context

 92:    Notes:
 93:    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.

 95:    If dynamic libraries are used, then the fourth input argument (routine_create)
 96:    is ignored.

 98:    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
 99:    and others of the form ${any_environmental_variable} occuring in pathname will be 
100:    replaced with appropriate values.

102:    Sample usage:
103: .vb
104:    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
105:                 "MySolverCreate",MySolverCreate);
106: .ve

108:    Then, your solver can be chosen with the procedural interface via
109: $     SNESSetType(snes,"my_solver")
110:    or at runtime via the option
111: $     -snes_type my_solver

113:    Level: advanced

115:     Note: If your function is not being put into a shared library then use SNESRegister() instead

117: .keywords: SNES, nonlinear, register

119: .seealso: SNESRegisterAll(), SNESRegisterDestroy()
120: M*/
121: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
122: #define SNESRegisterDynamic(a,b,c,d) SNESRegister(a,b,c,0)
123: #else
124: #define SNESRegisterDynamic(a,b,c,d) SNESRegister(a,b,c,d)
125: #endif

127: extern PetscErrorCode  SNESGetKSP(SNES,KSP*);
128: extern PetscErrorCode  SNESSetKSP(SNES,KSP);
129: extern PetscErrorCode  SNESGetSolution(SNES,Vec*);
130: extern PetscErrorCode  SNESGetSolutionUpdate(SNES,Vec*);
131: extern PetscErrorCode  SNESGetRhs(SNES,Vec*);
132: extern PetscErrorCode  SNESView(SNES,PetscViewer);

134: extern PetscErrorCode  SNESSetOptionsPrefix(SNES,const char[]);
135: extern PetscErrorCode  SNESAppendOptionsPrefix(SNES,const char[]);
136: extern PetscErrorCode  SNESGetOptionsPrefix(SNES,const char*[]);
137: extern PetscErrorCode  SNESSetFromOptions(SNES);
138: extern PetscErrorCode  SNESDefaultGetWork(SNES,PetscInt);

140: extern PetscErrorCode  MatCreateSNESMF(SNES,Mat*);
141: extern PetscErrorCode  MatMFFDComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);

143: extern PetscErrorCode  MatDAADSetSNES(Mat,SNES);

145: extern PetscErrorCode  SNESGetType(SNES,const SNESType*);
146: extern PetscErrorCode  SNESMonitorDefault(SNES,PetscInt,PetscReal,void *);
147: extern PetscErrorCode  SNESMonitorRange(SNES,PetscInt,PetscReal,void *);
148: extern PetscErrorCode  SNESMonitorRatio(SNES,PetscInt,PetscReal,void *);
149: extern PetscErrorCode  SNESMonitorSetRatio(SNES,PetscViewer);
150: extern PetscErrorCode  SNESMonitorSolution(SNES,PetscInt,PetscReal,void *);
151: extern PetscErrorCode  SNESMonitorResidual(SNES,PetscInt,PetscReal,void *);
152: extern PetscErrorCode  SNESMonitorSolutionUpdate(SNES,PetscInt,PetscReal,void *);
153: extern PetscErrorCode  SNESMonitorDefaultShort(SNES,PetscInt,PetscReal,void *);
154: extern PetscErrorCode  SNESSetTolerances(SNES,PetscReal,PetscReal,PetscReal,PetscInt,PetscInt);
155: extern PetscErrorCode  SNESGetTolerances(SNES,PetscReal*,PetscReal*,PetscReal*,PetscInt*,PetscInt*);
156: extern PetscErrorCode  SNESSetTrustRegionTolerance(SNES,PetscReal);
157: extern PetscErrorCode  SNESGetFunctionNorm(SNES,PetscReal*);
158: extern PetscErrorCode  SNESSetFunctionNorm(SNES,PetscReal);
159: extern PetscErrorCode  SNESGetIterationNumber(SNES,PetscInt*);
160: extern PetscErrorCode  SNESSetIterationNumber(SNES,PetscInt);

162: extern PetscErrorCode  SNESGetNonlinearStepFailures(SNES,PetscInt*);
163: extern PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES,PetscInt);
164: extern PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES,PetscInt*);
165: extern PetscErrorCode  SNESGetNumberFunctionEvals(SNES,PetscInt*);

167: extern PetscErrorCode  SNESSetLagPreconditioner(SNES,PetscInt);
168: extern PetscErrorCode  SNESGetLagPreconditioner(SNES,PetscInt*);
169: extern PetscErrorCode  SNESSetLagJacobian(SNES,PetscInt);
170: extern PetscErrorCode  SNESGetLagJacobian(SNES,PetscInt*);
171: extern PetscErrorCode  SNESSetGridSequence(SNES,PetscInt);

173: extern PetscErrorCode  SNESGetLinearSolveIterations(SNES,PetscInt*);
174: extern PetscErrorCode  SNESGetLinearSolveFailures(SNES,PetscInt*);
175: extern PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES,PetscInt);
176: extern PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES,PetscInt*);

178: extern PetscErrorCode  SNESKSPSetUseEW(SNES,PetscBool );
179: extern PetscErrorCode  SNESKSPGetUseEW(SNES,PetscBool *);
180: extern PetscErrorCode  SNESKSPSetParametersEW(SNES,PetscInt,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal);
181: extern PetscErrorCode  SNESKSPGetParametersEW(SNES,PetscInt*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*,PetscReal*);

183: extern PetscErrorCode  SNESMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
184: extern PetscErrorCode  SNESMonitorLG(SNES,PetscInt,PetscReal,void*);
185: extern PetscErrorCode  SNESMonitorLGDestroy(PetscDrawLG*);
186: extern PetscErrorCode  SNESMonitorLGRangeCreate(const char[],const char[],int,int,int,int,PetscDrawLG*);
187: extern PetscErrorCode  SNESMonitorLGRange(SNES,PetscInt,PetscReal,void*);
188: extern PetscErrorCode  SNESMonitorLGRangeDestroy(PetscDrawLG*);

190: extern PetscErrorCode  SNESSetApplicationContext(SNES,void *);
191: extern PetscErrorCode  SNESGetApplicationContext(SNES,void *);
192: extern PetscErrorCode  SNESSetComputeApplicationContext(SNES,PetscErrorCode (*)(SNES,void**),PetscErrorCode (*)(void**));

194: extern PetscErrorCode  SNESPythonSetType(SNES,const char[]);

196: extern PetscErrorCode  SNESSetFunctionDomainError(SNES);
197: /*E
198:     SNESConvergedReason - reason a SNES method was said to 
199:          have converged or diverged

201:    Level: beginner

203:    The two most common reasons for divergence are 
204: $   1) an incorrectly coded or computed Jacobian or 
205: $   2) failure or lack of convergence in the linear system (in this case we recommend
206: $      testing with -pc_type lu to eliminate the linear solver as the cause of the problem).

208:    Diverged Reasons:
209: .    SNES_DIVERGED_LOCAL_MIN - this can only occur when using the line-search variant of SNES.
210:        The line search wants to minimize Q(alpha) = 1/2 || F(x + alpha s) ||^2_2  this occurs
211:        at Q'(alpha) = s^T F'(x+alpha s)^T F(x+alpha s) = 0. If s is the Newton direction - F'(x)^(-1)F(x) then
212:        you get Q'(alpha) = -F(x)^T F'(x)^(-1)^T F'(x+alpha s)F(x+alpha s); when alpha = 0
213:        Q'(0) = - ||F(x)||^2_2 which is always NEGATIVE if F'(x) is invertible. This means the Newton
214:        direction is a descent direction and the line search should succeed if alpha is small enough.

216:        If F'(x) is NOT invertible AND F'(x)^T F(x) = 0 then Q'(0) = 0 and the Newton direction 
217:        is NOT a descent direction so the line search will fail. All one can do at this point
218:        is change the initial guess and try again.

220:        An alternative explanation: Newton's method can be regarded as replacing the function with
221:        its linear approximation and minimizing the 2-norm of that. That is F(x+s) approx F(x) + F'(x)s
222:        so we minimize || F(x) + F'(x) s ||^2_2; do this using Least Squares. If F'(x) is invertible then
223:        s = - F'(x)^(-1)F(x) otherwise F'(x)^T F'(x) s = -F'(x)^T F(x). If F'(x)^T F(x) is NOT zero then there
224:        exists a nontrival (that is F'(x)s != 0) solution to the equation and this direction is 
225:        s = - [F'(x)^T F'(x)]^(-1) F'(x)^T F(x) so Q'(0) = - F(x)^T F'(x) [F'(x)^T F'(x)]^(-T) F'(x)^T F(x)
226:        = - (F'(x)^T F(x)) [F'(x)^T F'(x)]^(-T) (F'(x)^T F(x)). Since we are assuming (F'(x)^T F(x)) != 0
227:        and F'(x)^T F'(x) has no negative eigenvalues Q'(0) < 0 so s is a descent direction and the line
228:        search should succeed for small enough alpha.

230:        Note that this RARELY happens in practice. Far more likely the linear system is not being solved
231:        (well enough?) or the Jacobian is wrong.
232:      
233:    SNES_DIVERGED_MAX_IT means that the solver reached the maximum number of iterations without satisfying any
234:    convergence criteria. SNES_CONVERGED_ITS means that SNESSkipConverged() was chosen as the convergence test;
235:    thus the usual convergence criteria have not been checked and may or may not be satisfied.

237:    Developer Notes: this must match finclude/petscsnes.h 

239:        The string versions of these are in SNESConvergedReason, if you change any value here you must
240:      also adjust that array.

242:    Each reason has its own manual page.

244: .seealso: SNESSolve(), SNESGetConvergedReason(), KSPConvergedReason, SNESSetConvergenceTest()
245: E*/
246: typedef enum {/* converged */
247:               SNES_CONVERGED_FNORM_ABS         =  2, /* ||F|| < atol */
248:               SNES_CONVERGED_FNORM_RELATIVE    =  3, /* ||F|| < rtol*||F_initial|| */
249:               SNES_CONVERGED_PNORM_RELATIVE    =  4, /* Newton computed step size small; || delta x || < stol */
250:               SNES_CONVERGED_ITS               =  5, /* maximum iterations reached */
251:               SNES_CONVERGED_TR_DELTA          =  7,
252:               /* diverged */
253:               SNES_DIVERGED_FUNCTION_DOMAIN     = -1, /* the new x location passed the function is not in the domain of F */
254:               SNES_DIVERGED_FUNCTION_COUNT      = -2,
255:               SNES_DIVERGED_LINEAR_SOLVE        = -3, /* the linear solve failed */
256:               SNES_DIVERGED_FNORM_NAN           = -4,
257:               SNES_DIVERGED_MAX_IT              = -5,
258:               SNES_DIVERGED_LINE_SEARCH         = -6, /* the line search failed */
259:               SNES_DIVERGED_INNER               = -7, /* inner solve failed */
260:               SNES_DIVERGED_LOCAL_MIN           = -8, /* || J^T b || is small, implies converged to local minimum of F() */
261:               SNES_CONVERGED_ITERATING          =  0} SNESConvergedReason;
262: extern const char *const*SNESConvergedReasons;

264: /*MC
265:      SNES_CONVERGED_FNORM_ABS - 2-norm(F) <= abstol

267:    Level: beginner

269: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

271: M*/

273: /*MC
274:      SNES_CONVERGED_FNORM_RELATIVE - 2-norm(F) <= rtol*2-norm(F(x_0)) where x_0 is the initial guess

276:    Level: beginner

278: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

280: M*/

282: /*MC
283:      SNES_CONVERGED_PNORM_RELATIVE - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
284:           solution and stol is the 4th argument to SNESSetTolerances()

286:    Level: beginner

288: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

290: M*/

292: /*MC
293:      SNES_DIVERGED_FUNCTION_COUNT - The user provided function has been called more times then the final
294:          argument to SNESSetTolerances()

296:    Level: beginner

298: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

300: M*/

302: /*MC
303:      SNES_DIVERGED_FNORM_NAN - the 2-norm of the current function evaluation is not-a-number (NaN), this
304:       is usually caused by a division of 0 by 0.

306:    Level: beginner

308: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

310: M*/

312: /*MC
313:      SNES_DIVERGED_MAX_IT - SNESSolve() has reached the maximum number of iterations requested

315:    Level: beginner

317: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

319: M*/

321: /*MC
322:      SNES_DIVERGED_LINE_SEARCH - The line search has failed. This only occurs for a SNESType of SNESLS

324:    Level: beginner

326: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

328: M*/

330: /*MC
331:      SNES_DIVERGED_LOCAL_MIN - the algorithm seems to have stagnated at a local minimum that is not zero. 
332:         See the manual page for SNESConvergedReason for more details

334:    Level: beginner

336: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

338: M*/

340: /*MC
341:      SNES_CONERGED_ITERATING - this only occurs if SNESGetConvergedReason() is called during the SNESSolve()

343:    Level: beginner

345: .seealso:  SNESSolve(), SNESGetConvergedReason(), SNESConvergedReason, SNESSetTolerances()

347: M*/

349: extern PetscErrorCode  SNESSetConvergenceTest(SNES,PetscErrorCode (*)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void*,PetscErrorCode (*)(void*));
350: extern PetscErrorCode  SNESDefaultConverged(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*);
351: extern PetscErrorCode  SNESSkipConverged(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*);
352: extern PetscErrorCode  SNESGetConvergedReason(SNES,SNESConvergedReason*);

354: extern PetscErrorCode  SNESDMDAComputeFunction(SNES,Vec,Vec,void*);
355: extern PetscErrorCode  SNESDMDAComputeJacobianWithAdic(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
356: extern PetscErrorCode  SNESDMDAComputeJacobianWithAdifor(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
357: extern PetscErrorCode  SNESDMDAComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);

359: extern PetscErrorCode SNESDMMeshComputeFunction(SNES,Vec,Vec,void*);
360: extern PetscErrorCode SNESDMMeshComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
361: extern PetscErrorCode SNESDMComplexComputeFunction(SNES,Vec,Vec,void *);
362: extern PetscErrorCode SNESDMComplexComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);

364: /* --------- Solving systems of nonlinear equations --------------- */
365: typedef PetscErrorCode (*SNESFunction)(SNES,Vec,Vec,void*);
366: typedef PetscErrorCode (*SNESJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
367: typedef PetscErrorCode (*SNESGSFunction)(SNES,Vec,Vec,void*);
368: extern PetscErrorCode  SNESSetFunction(SNES,Vec,SNESFunction,void*);
369: extern PetscErrorCode  SNESGetFunction(SNES,Vec*,SNESFunction*,void**);
370: extern PetscErrorCode  SNESComputeFunction(SNES,Vec,Vec);
371: extern PetscErrorCode  SNESSetJacobian(SNES,Mat,Mat,SNESJacobian,void*);
372: extern PetscErrorCode  SNESGetJacobian(SNES,Mat*,Mat*,SNESJacobian*,void**);
373: extern PetscErrorCode  SNESDefaultComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
374: extern PetscErrorCode  SNESDefaultComputeJacobianColor(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
375: extern PetscErrorCode  SNESSetComputeInitialGuess(SNES,PetscErrorCode (*)(SNES,Vec,void*),void*);
376: extern PetscErrorCode  SNESSetPicard(SNES,Vec,SNESFunction,Mat,Mat,SNESJacobian,void*);
377: extern PetscErrorCode  SNESGetPicard(SNES,Vec*,SNESFunction*,Mat*,SNESJacobian*,void**);

379: extern PetscErrorCode  SNESSetGS(SNES,SNESGSFunction,void*);
380: extern PetscErrorCode  SNESGetGS(SNES,SNESGSFunction*,void**);
381: extern PetscErrorCode  SNESSetUseGS(SNES,PetscBool);
382: extern PetscErrorCode  SNESGetUseGS(SNES,PetscBool *);
383: extern PetscErrorCode  SNESSetGSSweeps(SNES,PetscInt);
384: extern PetscErrorCode  SNESGetGSSweeps(SNES,PetscInt *);
385: extern PetscErrorCode  SNESComputeGS(SNES,Vec,Vec);

387: /* --------- Routines specifically for line search methods --------------- */
388: /*E
389:     SNESLineSearchType - type of line search used in Newton's method as well as VI solvers and Richardson solvers

391:     Level: beginner

393: .seealso: SNESSetFromOptions(), SNESLineSearchSet()
394: E*/
395: typedef enum {SNES_LS_BASIC, SNES_LS_BASIC_NONORMS, SNES_LS_QUADRATIC, SNES_LS_CUBIC, SNES_LS_EXACT, SNES_LS_TEST, SNES_LS_SECANT,SNES_LS_USER_DEFINED} SNESLineSearchType;
396: extern const char *const SNESLineSearchTypes[];
397: extern const char *SNESLineSearchTypeName(SNESLineSearchType); /* Does bounds checking, use this for viewing */

399: extern PetscErrorCode  SNESLineSearchSet(SNES,PetscErrorCode(*)(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *),void*);
400: extern PetscErrorCode  SNESLineSearchSetType(SNES,SNESLineSearchType);
401: extern PetscErrorCode  SNESLineSearchNo(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *);
402: extern PetscErrorCode  SNESLineSearchNoNorms(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *);
403: extern PetscErrorCode  SNESLineSearchQuadratic(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *);
404: extern PetscErrorCode  SNESLineSearchCubic(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *);
405: extern PetscErrorCode  SNESLineSearchSecant(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *);
406: extern PetscErrorCode  SNESLineSearchQuadraticSecant(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *);

408: extern PetscErrorCode  SNESLineSearchApply(SNES,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *);
409: extern PetscErrorCode  SNESLineSearchPreCheckApply(SNES,Vec,Vec,PetscBool*);
410: extern PetscErrorCode  SNESLineSearchPostCheckApply(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*);

412: extern PetscErrorCode  SNESLineSearchSetPostCheck(SNES,PetscErrorCode(*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*);
413: extern PetscErrorCode  SNESLineSearchSetPreCheck(SNES,PetscErrorCode(*)(SNES,Vec,Vec,void*,PetscBool *),void*);
414: extern PetscErrorCode  SNESLineSearchPreCheckPicard(SNES,Vec,Vec,void*,PetscBool*);
415: extern PetscErrorCode  SNESLineSearchSetParams(SNES,PetscReal,PetscReal,PetscReal);
416: extern PetscErrorCode  SNESLineSearchGetParams(SNES,PetscReal*,PetscReal*,PetscReal*);
417: extern PetscErrorCode  SNESLineSearchSetMonitor(SNES,PetscBool );

419: extern PetscErrorCode  SNESShellGetContext(SNES,void**);
420: extern PetscErrorCode  SNESShellSetContext(SNES,void*);
421: extern PetscErrorCode  SNESShellSetSolve(SNES,PetscErrorCode (*)(SNES,Vec));

423: /* Routines for VI solver */
424: extern PetscErrorCode  SNESVISetVariableBounds(SNES,Vec,Vec);
425: extern PetscErrorCode  SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES,Vec,Vec));
426: extern PetscErrorCode  SNESVIGetInactiveSet(SNES,IS*);
427: extern PetscErrorCode  SNESVIGetActiveSetIS(SNES,Vec,Vec,IS*);
428: extern PetscErrorCode  SNESVIComputeInactiveSetFnorm(SNES,Vec,Vec,PetscReal*);
429: extern PetscErrorCode  SNESVISetRedundancyCheck(SNES,PetscErrorCode(*)(SNES,IS,IS*,void*),void*);
430: #define SNES_VI_INF   1.0e20
431: #define SNES_VI_NINF -1.0e20

433: extern PetscErrorCode  SNESTestLocalMin(SNES);

435: /* Should this routine be private? */
436: extern PetscErrorCode  SNESComputeJacobian(SNES,Vec,Mat*,Mat*,MatStructure*);

438: extern PetscErrorCode SNESSetDM(SNES,DM);
439: extern PetscErrorCode SNESGetDM(SNES,DM*);
440: extern PetscErrorCode SNESSetPC(SNES,SNES);
441: extern PetscErrorCode SNESGetPC(SNES,SNES*);
442: extern PetscErrorCode SNESRestrictHookAdd(SNES,PetscErrorCode (*)(SNES,SNES,void*),void*);
443: extern PetscErrorCode SNESRestrictHooksRun(SNES,SNES);

445: /* Routines for Multiblock solver */
446: extern PetscErrorCode SNESMultiblockSetFields(SNES, const char [], PetscInt, const PetscInt *);
447: extern PetscErrorCode SNESMultiblockSetIS(SNES, const char [], IS);
448: extern PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
449: extern PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);

451: /*J
452:     SNESMSType - String with the name of a PETSc SNESMS method.

454:    Level: intermediate

456: .seealso: SNESMSSetType(), SNES
457: J*/
458: #define SNESMSType char*
459: #define SNESMSM62       "m62"
460: #define SNESMSEULER     "euler"
461: #define SNESMSJAMESON83 "jameson83"
462: #define SNESMSVLTP21    "vltp21"
463: #define SNESMSVLTP31    "vltp31"
464: #define SNESMSVLTP41    "vltp41"
465: #define SNESMSVLTP51    "vltp51"
466: #define SNESMSVLTP61    "vltp61"

468: extern PetscErrorCode SNESMSRegister(const SNESMSType,PetscInt,PetscInt,PetscReal,const PetscReal[],const PetscReal[],const PetscReal[]);
469: extern PetscErrorCode SNESMSSetType(SNES,const SNESMSType);
470: extern PetscErrorCode SNESMSFinalizePackage(void);
471: extern PetscErrorCode SNESMSInitializePackage(const char path[]);
472: extern PetscErrorCode SNESMSRegisterDestroy(void);
473: extern PetscErrorCode SNESMSRegisterAll(void);

475: PETSC_EXTERN_CXX_END
476: #endif