Actual source code: petscsnes.h

  1: /*
  2:     User interface for the nonlinear solvers package.
  3: */
  4: #pragma once

  6: #include <petscksp.h>
  7: #include <petscdmtypes.h>
  8: #include <petscfvtypes.h>
  9: #include <petscdmdatypes.h>
 10: #include <petscsnestypes.h>

 12: /* SUBMANSEC = SNES */

 14: /*J
 15:    SNESType - String with the name of a PETSc `SNES` method. These are all the nonlinear solvers that PETSc provides.

 17:    Level: beginner

 19:    Note:
 20:    Use `SNESSetType()` or the options database key `-snes_type` to set the specific nonlinear solver algorithm to use with a given `SNES` object

 22: .seealso: [](doc_nonlinsolve), [](ch_snes), `SNESSetType()`, `SNES`, `SNESCreate()`, `SNESDestroy()`, `SNESSetFromOptions()`
 23: J*/
 24: typedef const char *SNESType;
 25: #define SNESNEWTONLS         "newtonls"
 26: #define SNESNEWTONTR         "newtontr"
 27: #define SNESNEWTONTRDC       "newtontrdc"
 28: #define SNESPYTHON           "python"
 29: #define SNESNRICHARDSON      "nrichardson"
 30: #define SNESKSPONLY          "ksponly"
 31: #define SNESKSPTRANSPOSEONLY "ksptransposeonly"
 32: #define SNESVINEWTONRSLS     "vinewtonrsls"
 33: #define SNESVINEWTONSSLS     "vinewtonssls"
 34: #define SNESNGMRES           "ngmres"
 35: #define SNESQN               "qn"
 36: #define SNESSHELL            "shell"
 37: #define SNESNGS              "ngs"
 38: #define SNESNCG              "ncg"
 39: #define SNESFAS              "fas"
 40: #define SNESMS               "ms"
 41: #define SNESNASM             "nasm"
 42: #define SNESANDERSON         "anderson"
 43: #define SNESASPIN            "aspin"
 44: #define SNESCOMPOSITE        "composite"
 45: #define SNESPATCH            "patch"
 46: #define SNESNEWTONAL         "newtonal"

 48: /* Logging support */
 49: PETSC_EXTERN PetscClassId SNES_CLASSID;
 50: PETSC_EXTERN PetscClassId DMSNES_CLASSID;

 52: PETSC_EXTERN PetscErrorCode SNESInitializePackage(void);
 53: PETSC_EXTERN PetscErrorCode SNESFinalizePackage(void);

 55: PETSC_EXTERN PetscErrorCode SNESCreate(MPI_Comm, SNES *);
 56: PETSC_EXTERN PetscErrorCode SNESParametersInitialize(SNES);
 57: PETSC_EXTERN PetscErrorCode SNESReset(SNES);
 58: PETSC_EXTERN PetscErrorCode SNESDestroy(SNES *);
 59: PETSC_EXTERN PetscErrorCode SNESSetType(SNES, SNESType);
 60: PETSC_EXTERN PetscErrorCode SNESMonitor(SNES, PetscInt, PetscReal);
 61: PETSC_EXTERN PetscErrorCode SNESMonitorSet(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *), void *, PetscCtxDestroyFn *);
 62: PETSC_EXTERN PetscErrorCode SNESMonitorSetFromOptions(SNES, const char[], const char[], const char[], PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(SNES, PetscViewerAndFormat *));
 63: PETSC_EXTERN PetscErrorCode SNESMonitorCancel(SNES);
 64: PETSC_EXTERN PetscErrorCode SNESMonitorSAWs(SNES, PetscInt, PetscReal, void *);
 65: PETSC_EXTERN PetscErrorCode SNESMonitorSAWsCreate(SNES, void **);
 66: PETSC_EXTERN PetscErrorCode SNESMonitorSAWsDestroy(void **);
 67: PETSC_EXTERN PetscErrorCode SNESSetConvergenceHistory(SNES, PetscReal[], PetscInt[], PetscInt, PetscBool);
 68: PETSC_EXTERN PetscErrorCode SNESGetConvergenceHistory(SNES, PetscReal *[], PetscInt *[], PetscInt *);
 69: PETSC_EXTERN PetscErrorCode SNESSetUp(SNES);
 70: PETSC_EXTERN PetscErrorCode SNESSolve(SNES, Vec, Vec);
 71: PETSC_EXTERN PetscErrorCode SNESSetErrorIfNotConverged(SNES, PetscBool);
 72: PETSC_EXTERN PetscErrorCode SNESGetErrorIfNotConverged(SNES, PetscBool *);
 73: PETSC_EXTERN PetscErrorCode SNESConverged(SNES, PetscInt, PetscReal, PetscReal, PetscReal);

 75: PETSC_EXTERN PetscErrorCode SNESSetWorkVecs(SNES, PetscInt);

 77: PETSC_EXTERN PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*)(SNES));

 79: PETSC_EXTERN PetscErrorCode SNESRegister(const char[], PetscErrorCode (*)(SNES));

 81: PETSC_EXTERN PetscErrorCode SNESGetKSP(SNES, KSP *);
 82: PETSC_EXTERN PetscErrorCode SNESSetKSP(SNES, KSP);
 83: PETSC_EXTERN PetscErrorCode SNESSetSolution(SNES, Vec);
 84: PETSC_EXTERN PetscErrorCode SNESGetSolution(SNES, Vec *);
 85: PETSC_EXTERN PetscErrorCode SNESGetSolutionUpdate(SNES, Vec *);
 86: PETSC_EXTERN PetscErrorCode SNESGetRhs(SNES, Vec *);
 87: PETSC_EXTERN PetscErrorCode SNESView(SNES, PetscViewer);
 88: PETSC_EXTERN PetscErrorCode SNESLoad(SNES, PetscViewer);
 89: PETSC_EXTERN PetscErrorCode SNESConvergedReasonViewSet(SNES, PetscErrorCode (*)(SNES, void *), void *, PetscCtxDestroyFn *);
 90: PETSC_EXTERN PetscErrorCode SNESViewFromOptions(SNES, PetscObject, const char[]);
 91: PETSC_EXTERN PetscErrorCode SNESConvergedReasonView(SNES, PetscViewer);
 92: PETSC_EXTERN PetscErrorCode SNESConvergedReasonViewFromOptions(SNES);
 93: PETSC_EXTERN PetscErrorCode SNESConvergedReasonViewCancel(SNES);

 95: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SNESConvergedReasonView()", ) static inline PetscErrorCode SNESReasonView(SNES snes, PetscViewer v)
 96: {
 97:   return SNESConvergedReasonView(snes, v);
 98: }
 99: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SNESConvergedReasonViewFromOptions()", ) static inline PetscErrorCode SNESReasonViewFromOptions(SNES snes)
100: {
101:   return SNESConvergedReasonViewFromOptions(snes);
102: }

104: #define SNES_FILE_CLASSID 1211224

106: PETSC_EXTERN PetscErrorCode SNESSetOptionsPrefix(SNES, const char[]);
107: PETSC_EXTERN PetscErrorCode SNESAppendOptionsPrefix(SNES, const char[]);
108: PETSC_EXTERN PetscErrorCode SNESGetOptionsPrefix(SNES, const char *[]);
109: PETSC_EXTERN PetscErrorCode SNESSetFromOptions(SNES);
110: PETSC_EXTERN PetscErrorCode SNESResetFromOptions(SNES);

112: PETSC_EXTERN PetscErrorCode SNESSetUseMatrixFree(SNES, PetscBool, PetscBool);
113: PETSC_EXTERN PetscErrorCode SNESGetUseMatrixFree(SNES, PetscBool *, PetscBool *);
114: PETSC_EXTERN PetscErrorCode MatCreateSNESMF(SNES, Mat *);
115: PETSC_EXTERN PetscErrorCode MatSNESMFGetSNES(Mat, SNES *);
116: PETSC_EXTERN PetscErrorCode MatSNESMFSetReuseBase(Mat, PetscBool);
117: PETSC_EXTERN PetscErrorCode MatSNESMFGetReuseBase(Mat, PetscBool *);
118: PETSC_EXTERN PetscErrorCode MatMFFDComputeJacobian(SNES, Vec, Mat, Mat, void *);
119: PETSC_EXTERN PetscErrorCode MatCreateSNESMFMore(SNES, Vec, Mat *);
120: PETSC_EXTERN PetscErrorCode MatSNESMFMoreSetParameters(Mat, PetscReal, PetscReal, PetscReal);

122: PETSC_EXTERN PetscErrorCode SNESGetType(SNES, SNESType *);
123: PETSC_EXTERN PetscErrorCode SNESMonitorDefaultSetUp(SNES, PetscViewerAndFormat *);
124: PETSC_EXTERN PetscErrorCode SNESMonitorDefault(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
125: PETSC_EXTERN PetscErrorCode SNESMonitorScaling(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
126: PETSC_EXTERN PetscErrorCode SNESMonitorRange(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
127: PETSC_EXTERN PetscErrorCode SNESMonitorRatio(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
128: PETSC_EXTERN PetscErrorCode SNESMonitorRatioSetUp(SNES, PetscViewerAndFormat *);
129: PETSC_EXTERN PetscErrorCode SNESMonitorSolution(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
130: PETSC_EXTERN PetscErrorCode SNESMonitorResidual(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
131: PETSC_EXTERN PetscErrorCode SNESMonitorSolutionUpdate(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
132: PETSC_EXTERN PetscErrorCode SNESMonitorDefaultShort(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
133: PETSC_EXTERN PetscErrorCode SNESMonitorDefaultField(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
134: PETSC_EXTERN PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
135: PETSC_EXTERN PetscErrorCode SNESMonitorFields(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
136: PETSC_EXTERN PetscErrorCode KSPMonitorSNESResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
137: PETSC_EXTERN PetscErrorCode KSPMonitorSNESResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
138: PETSC_EXTERN PetscErrorCode KSPMonitorSNESResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);

140: PETSC_EXTERN PetscErrorCode SNESSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt, PetscInt);
141: PETSC_EXTERN PetscErrorCode SNESSetDivergenceTolerance(SNES, PetscReal);
142: PETSC_EXTERN PetscErrorCode SNESGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *, PetscInt *);
143: PETSC_EXTERN PetscErrorCode SNESGetDivergenceTolerance(SNES, PetscReal *);
144: PETSC_EXTERN PetscErrorCode SNESGetForceIteration(SNES, PetscBool *);
145: PETSC_EXTERN PetscErrorCode SNESSetForceIteration(SNES, PetscBool);
146: PETSC_EXTERN PetscErrorCode SNESGetIterationNumber(SNES, PetscInt *);
147: PETSC_EXTERN PetscErrorCode SNESSetIterationNumber(SNES, PetscInt);

149: /*E
150:    SNESNewtonTRFallbackType - type of fallback in case the solution of the trust-region subproblem is outside of the radius

152:    Values:
153: +  `SNES_TR_FALLBACK_NEWTON` - use scaled Newton step
154: .  `SNES_TR_FALLBACK_CAUCHY` - use Cauchy direction
155: -  `SNES_TR_FALLBACK_DOGLEG` - use dogleg method

157:    Level: intermediate

159: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNEWTONTRDC`
160: E*/
161: typedef enum {
162:   SNES_TR_FALLBACK_NEWTON,
163:   SNES_TR_FALLBACK_CAUCHY,
164:   SNES_TR_FALLBACK_DOGLEG,
165: } SNESNewtonTRFallbackType;

167: PETSC_EXTERN const char *const SNESNewtonTRFallbackTypes[];

169: PETSC_EXTERN PetscErrorCode SNESNewtonTRSetPreCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, PetscBool *, void *), void *ctx);
170: PETSC_EXTERN PetscErrorCode SNESNewtonTRGetPreCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, PetscBool *, void *), void **ctx);
171: PETSC_EXTERN PetscErrorCode SNESNewtonTRSetPostCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);
172: PETSC_EXTERN PetscErrorCode SNESNewtonTRGetPostCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);
173: PETSC_EXTERN PetscErrorCode SNESNewtonTRSetFallbackType(SNES, SNESNewtonTRFallbackType);
174: PETSC_EXTERN PetscErrorCode SNESNewtonTRPreCheck(SNES, Vec, Vec, PetscBool *);
175: PETSC_EXTERN PetscErrorCode SNESNewtonTRPostCheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *);
176: PETSC_EXTERN PetscErrorCode SNESNewtonTRSetNormType(SNES, NormType);

178: /*E
179:     SNESNewtonTRQNType - type of quasi-Newton model to use

181:    Values:
182: +  `SNES_TR_QN_NONE` - do not use a quasi-Newton model
183: .  `SNES_TR_QN_SAME` - use the same quasi-Newton model for matrix and preconditioner
184: -  `SNES_TR_QN_DIFFERENT` - use different quasi-Newton models for matrix and preconditioner

186:    Level: intermediate

188: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`
189: E*/
190: typedef enum {
191:   SNES_TR_QN_NONE,
192:   SNES_TR_QN_SAME,
193:   SNES_TR_QN_DIFFERENT,
194: } SNESNewtonTRQNType;

196: PETSC_EXTERN const char *const SNESNewtonTRQNTypes[];

198: PETSC_EXTERN PetscErrorCode SNESNewtonTRSetQNType(SNES, SNESNewtonTRQNType);

200: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 22, 0, "SNESNewtonTRSetTolerances()", ) PetscErrorCode SNESSetTrustRegionTolerance(SNES, PetscReal);
201: PETSC_EXTERN PetscErrorCode SNESNewtonTRSetTolerances(SNES, PetscReal, PetscReal, PetscReal);
202: PETSC_EXTERN PetscErrorCode SNESNewtonTRGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *);
203: PETSC_EXTERN PetscErrorCode SNESNewtonTRSetUpdateParameters(SNES, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
204: PETSC_EXTERN PetscErrorCode SNESNewtonTRGetUpdateParameters(SNES, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *);

206: PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES, PetscBool *);
207: PETSC_EXTERN PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, PetscBool *, void *), void *ctx);
208: PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, PetscBool *, void *), void **ctx);
209: PETSC_EXTERN PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);
210: PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);

212: PETSC_EXTERN PetscErrorCode SNESGetNonlinearStepFailures(SNES, PetscInt *);
213: PETSC_EXTERN PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES, PetscInt);
214: PETSC_EXTERN PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES, PetscInt *);
215: PETSC_EXTERN PetscErrorCode SNESGetNumberFunctionEvals(SNES, PetscInt *);

217: PETSC_EXTERN PetscErrorCode SNESSetLagPreconditioner(SNES, PetscInt);
218: PETSC_EXTERN PetscErrorCode SNESGetLagPreconditioner(SNES, PetscInt *);
219: PETSC_EXTERN PetscErrorCode SNESSetLagJacobian(SNES, PetscInt);
220: PETSC_EXTERN PetscErrorCode SNESGetLagJacobian(SNES, PetscInt *);
221: PETSC_EXTERN PetscErrorCode SNESSetLagPreconditionerPersists(SNES, PetscBool);
222: PETSC_EXTERN PetscErrorCode SNESSetLagJacobianPersists(SNES, PetscBool);
223: PETSC_EXTERN PetscErrorCode SNESSetGridSequence(SNES, PetscInt);
224: PETSC_EXTERN PetscErrorCode SNESGetGridSequence(SNES, PetscInt *);

226: PETSC_EXTERN PetscErrorCode SNESGetLinearSolveIterations(SNES, PetscInt *);
227: PETSC_EXTERN PetscErrorCode SNESGetLinearSolveFailures(SNES, PetscInt *);
228: PETSC_EXTERN PetscErrorCode SNESSetMaxLinearSolveFailures(SNES, PetscInt);
229: PETSC_EXTERN PetscErrorCode SNESGetMaxLinearSolveFailures(SNES, PetscInt *);
230: PETSC_EXTERN PetscErrorCode SNESSetCountersReset(SNES, PetscBool);
231: PETSC_EXTERN PetscErrorCode SNESResetCounters(SNES);

233: PETSC_EXTERN PetscErrorCode SNESKSPSetUseEW(SNES, PetscBool);
234: PETSC_EXTERN PetscErrorCode SNESKSPGetUseEW(SNES, PetscBool *);
235: PETSC_EXTERN PetscErrorCode SNESKSPSetParametersEW(SNES, PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
236: PETSC_EXTERN PetscErrorCode SNESKSPGetParametersEW(SNES, PetscInt *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *);

238: PETSC_EXTERN PetscErrorCode SNESMonitorLGRange(SNES, PetscInt, PetscReal, void *);

240: PETSC_EXTERN PetscErrorCode SNESSetApplicationContext(SNES, void *);
241: PETSC_EXTERN PetscErrorCode SNESGetApplicationContext(SNES, void *);
242: PETSC_EXTERN PetscErrorCode SNESSetComputeApplicationContext(SNES, PetscErrorCode (*)(SNES, void **), PetscCtxDestroyFn *);

244: PETSC_EXTERN PetscErrorCode SNESPythonSetType(SNES, const char[]);
245: PETSC_EXTERN PetscErrorCode SNESPythonGetType(SNES, const char *[]);

247: PETSC_EXTERN PetscErrorCode SNESSetFunctionDomainError(SNES);
248: PETSC_EXTERN PetscErrorCode SNESSetObjectiveDomainError(SNES);
249: PETSC_EXTERN PetscErrorCode SNESSetJacobianDomainError(SNES);
250: PETSC_EXTERN PetscErrorCode SNESSetCheckJacobianDomainError(SNES, PetscBool);
251: PETSC_EXTERN PetscErrorCode SNESGetCheckJacobianDomainError(SNES, PetscBool *);

253: #define SNES_CONVERGED_TR_DELTA_DEPRECATED SNES_CONVERGED_TR_DELTA PETSC_DEPRECATED_ENUM(3, 12, 0, "SNES_DIVERGED_TR_DELTA", )
254: /*E
255:     SNESConvergedReason - reason a `SNESSolve()` was determined to have converged or diverged

257:    Values:
258: +  `SNES_CONVERGED_FNORM_ABS`         - $ ||F|| \le abstol $
259: .  `SNES_CONVERGED_FNORM_RELATIVE`    - $ ||F|| <= rtol*||F(x_0))|| $ where $x_0 $ is the initial guess
260: .  `SNES_CONVERGED_SNORM_RELATIVE`    - The 2-norm of the last step $ \le stol * ||x|| $ where $ x $ is the current solution
261: .  `SNES_CONVERGED_USER`              - The user has indicated convergence for an arbitrary reason
262: .  `SNES_DIVERGED_FUNCTION_COUNT`     - The user provided function has been called more times than the maximum set in `SNESSetTolerances()`
263: .  `SNES_DIVERGED_DTOL`               - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
264: .  `SNES_DIVERGED_FUNCTION_NANORINF`  - the 2-norm of the current function evaluation is not-a-number (NaN) or infinity, (this
265:                                         is usually caused by a division of 0 by 0) and the solver could not recover from this (by, for example, cutting the step size)
266: .  `SNES_DIVERGED_OBJECTIVE_NANORINF` - the object function evaluation is not-a-number (NaN) or infinity, (this
267:                                         is usually caused by a division of 0 by 0) and the solver could not recover from this (by, for example, cutting the step size)
268: .  `SNES_DIVERGED_FUNCTION_DOMAIN`    - the function evaluation occurred outside the function's domain (function callback provided by
269:                                         `SNESSetFunction()` called `SNESSetObjectiveDomainError()`) and the solver could not recover from this (by, for example, cutting the step size)
270: .  `SNES_DIVERGED_OBJECTIVE_DOMAIN`   - the object function evaluation occurred outside the function's domain (function callback provided by
271:                                         `SNESSetObjective()` called `SNESSetObjectiveDomainError()`) and the solver could not recover from this (by, for example, cutting the step size)
272: .  `SNES_DIVERGED_JACOBIAN_DOMAIN`    - the Jacobian evaluation occurred outside the function's domain (function callback provided by
273:                                         `SNESSetJacobian()` called `SNESSetJacobianDomainError()`)
274: .  `SNES_DIVERGED_MAX_IT`             - `SNESSolve()` has reached the maximum number of iterations requested
275: .  `SNES_DIVERGED_LINE_SEARCH`        - The line search has failed. This only occurs for `SNES` solvers that use a line search
276: .  `SNES_DIVERGED_LOCAL_MIN`          - the algorithm seems to have stagnated at a local minimum that is not zero.
277: -  `SNES_CONVERGED_ITERATING          - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`

279:    Level: beginner

281:     Notes:
282:    The two most common reasons for divergence are an incorrectly coded or computed Jacobian or failure or lack of convergence in the linear system
283:    (in this case we recommend
284:    testing with `-pc_type lu` to eliminate the linear solver as the cause of the problem).

286:    `SNES_DIVERGED_LOCAL_MIN` can only occur when using a `SNES` solver that uses a line search (`SNESLineSearch`).
287:    The line search wants to $ \min Q(\alpha) = 1/2 || F(x + \alpha s) ||^2_2 $  this occurs
288:    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
289:    $ Q'(\alpha) = -F(x)^T F'(x)^(-1)^T F'(x+\alpha s)F(x+\alpha s)$; when $\alpha = 0$
290:    $Q'(0) = - ||F(x)||^2_2 $ which is always NEGATIVE if $F'(x)$ is invertible. This means the Newton
291:    direction is a descent direction and the line search should succeed if $\alpha $ is small enough.

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

297:    An alternative explanation: Newton's method can be regarded as replacing the function with
298:    its linear approximation and minimizing the 2-norm of that. That is $F(x+s) \approx F(x) + F'(x)s$
299:    so we minimize $ || F(x) + F'(x) s ||^2_2$ using Least Squares. If $F'(x)$ is invertible then
300:    $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
301:    exists a nontrivial (that is $F'(x)s \ne 0$) solution to the equation and this direction is
302:    $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)
303:    = - (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)) \ne 0$
304:    and $F'(x)^T F'(x)$ has no negative eigenvalues $Q'(0) < 0$ so $s$ is a descent direction and the line
305:    search should succeed for small enough $\alpha$.

307:    Note that this RARELY happens in practice. Far more likely the linear system is not being solved
308:    (well enough?) or the Jacobian is wrong.

310:    `SNES_DIVERGED_MAX_IT` means that the solver reached the maximum number of iterations without satisfying any
311:    convergence criteria. `SNES_CONVERGED_ITS` means that `SNESConvergedSkip()` was chosen as the convergence test;
312:    thus the usual convergence criteria have not been checked and may or may not be satisfied.

314: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`, `SNESSetTolerances()`
315: E*/
316: typedef enum {                       /* converged */
317:   SNES_CONVERGED_FNORM_ABS      = 2, /* ||F|| < atol */
318:   SNES_CONVERGED_FNORM_RELATIVE = 3, /* ||F|| < rtol*||F_initial|| */
319:   SNES_CONVERGED_SNORM_RELATIVE = 4, /* Newton computed step size small; || delta x || < stol || x || */
320:   SNES_CONVERGED_ITS            = 5, /* maximum iterations reached */
321:   SNES_BREAKOUT_INNER_ITER      = 6, /* Flag to break out of inner loop after checking custom convergence, used in multi-phase flow when state changes */
322:   SNES_CONVERGED_USER           = 7, /* The user has indicated convergence for an arbitrary reason */
323:   /* diverged */
324:   SNES_DIVERGED_FUNCTION_DOMAIN      = -1, /* the new x location passed the function is not in the domain of F */
325:   SNES_DIVERGED_FUNCTION_COUNT       = -2,
326:   SNES_DIVERGED_LINEAR_SOLVE         = -3, /* the linear solve failed */
327:   SNES_DIVERGED_FUNCTION_NANORINF    = -4,
328:   SNES_DIVERGED_MAX_IT               = -5,
329:   SNES_DIVERGED_LINE_SEARCH          = -6,  /* the line search failed */
330:   SNES_DIVERGED_INNER                = -7,  /* inner solve failed */
331:   SNES_DIVERGED_LOCAL_MIN            = -8,  /* || J^T b || is small, implies converged to local minimum of F() */
332:   SNES_DIVERGED_DTOL                 = -9,  /* || F || > divtol*||F_initial|| */
333:   SNES_DIVERGED_JACOBIAN_DOMAIN      = -10, /* Jacobian calculation does not make sense */
334:   SNES_DIVERGED_TR_DELTA             = -11,
335:   SNES_CONVERGED_TR_DELTA_DEPRECATED = -11,
336:   SNES_DIVERGED_USER                 = -12, /* The user has indicated divergence for an arbitrary reason */
337:   SNES_DIVERGED_OBJECTIVE_DOMAIN     = -13,
338:   SNES_DIVERGED_OBJECTIVE_NANORINF   = -14,

340:   SNES_CONVERGED_ITERATING = 0
341: } SNESConvergedReason;
342: PETSC_EXTERN const char *const *SNESConvergedReasons;

344: /*MC
345:    SNES_CONVERGED_FNORM_ABS - $||F|| \le abstol$

347:    Level: beginner

349: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
350: M*/

352: /*MC
353:    SNES_CONVERGED_FNORM_RELATIVE - $||F|| \le rtol*||F(x_0)||$ where $x_0$ is the initial guess

355:    Level: beginner

357: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
358: M*/

360: /*MC
361:   SNES_CONVERGED_SNORM_RELATIVE - The 2-norm of the last step $\le stol * ||x||$ where $x$ is the current
362:   solution and $stol$ is the 4th argument to `SNESSetTolerances()`

364:   Options Database Key:
365:   -snes_stol <stol> - the step tolerance

367:    Level: beginner

369: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
370: M*/

372: /*MC
373:    SNES_DIVERGED_FUNCTION_COUNT - The user provided function has been called more times then the final
374:    argument to `SNESSetTolerances()`

376:    Level: beginner

378: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
379: M*/

381: /*MC
382:    SNES_DIVERGED_DTOL - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`

384:    Level: beginner

386: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESSetDivergenceTolerance()`
387: M*/

389: /*MC
390:    SNES_DIVERGED_FUNCTION_NANORINF - the 2-norm of the current function evaluation is not-a-number (NaN) or infinity, this
391:    is usually caused by a division of 0 by 0, or infinity.  See `SNESSetFunctionDomainError()`

393:    Level: beginner

395: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
396: M*/

398: /*MC
399:    SNES_DIVERGED_FUNCTION_DOMAIN - the function provided with `SNESSetFunction()` called `SNESSetFunctionDomainError()` and
400:    the solver could not recoverer.

402:    Level: beginner

404: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
405: M*/

407: /*MC
408:    SNES_DIVERGED_OBJECTIVE_DOMAIN - the function provided with `SNESSetObjective()` called `SNESSetObjectiveDomainError()` and
409:    the solver could not recoverer.

411:    Level: beginner

413: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
414: M*/

416: /*MC
417:    SNES_DIVERGED_JACOBIAN_DOMAIN - the function provided with `SNESSetJacobian()` called `SNESSetJacobianDomainError()`

419:    Level: beginner

421: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
422: M*/

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

427:    Level: beginner

429: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
430: M*/

432: /*MC
433:    SNES_DIVERGED_LINE_SEARCH - The line search has failed. This only occurs for a `SNES` solvers that use a line search

435:    Level: beginner

437: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESLineSearch`
438: M*/

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

444:    Level: beginner

446: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
447: M*/

449: /*MC
450:    SNES_CONERGED_ITERATING - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`

452:    Level: beginner

454: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
455: M*/

457: PETSC_EXTERN PetscErrorCode SNESSetConvergenceTest(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *), void *, PetscCtxDestroyFn *);
458: PETSC_EXTERN PetscErrorCode SNESConvergedDefault(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
459: PETSC_EXTERN PetscErrorCode SNESConvergedSkip(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
460: PETSC_EXTERN PetscErrorCode SNESConvergedCorrectPressure(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
461: PETSC_EXTERN PetscErrorCode SNESGetConvergedReason(SNES, SNESConvergedReason *);
462: PETSC_EXTERN PetscErrorCode SNESGetConvergedReasonString(SNES, const char **);
463: PETSC_EXTERN PetscErrorCode SNESSetConvergedReason(SNES, SNESConvergedReason);

465: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "SNESConvergedSkip()", ) static inline void SNESSkipConverged(void)
466: { /* never called */
467: }
468: #define SNESSkipConverged (SNESSkipConverged, SNESConvergedSkip)

470: /*S
471:   SNESInitialGuessFn - A prototype of a `SNES` compute initial guess function that would be passed to `SNESSetComputeInitialGuess()`

473:   Calling Sequence:
474: + snes  - `SNES` context
475: . u   - output vector to contain initial guess
476: - ctx - [optional] user-defined function context

478:   Level: beginner

480: .seealso: [](ch_snes), `SNES`, `SNESSetComputeInitialGuess()`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESFunctionFn`
481: S*/
482: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESInitialGuessFn(SNES snes, Vec u, void *ctx);

484: /*S
485:   SNESFunctionFn - A prototype of a `SNES` evaluation function that would be passed to `SNESSetFunction()`

487:   Calling Sequence:
488: + snes  - `SNES` context
489: . u   - input vector
490: . F   - function vector
491: - ctx - [optional] user-defined function context

493:   Level: beginner

495: .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
496: S*/
497: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESFunctionFn(SNES snes, Vec u, Vec F, void *ctx);

499: /*S
500:   SNESObjectiveFn - A prototype of a `SNES` objective evaluation function that would be passed to `SNESSetObjective()`

502:   Calling Sequence:
503: + snes  - `SNES` context
504: . u   - input vector
505: . o   - output value
506: - ctx - [optional] user-defined function context

508:   Level: beginner

510: .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
511: S*/
512: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESObjectiveFn(SNES snes, Vec u, PetscReal *o, void *ctx);

514: /*S
515:   SNESJacobianFn - A prototype of a `SNES` Jacobian evaluation function that would be passed to `SNESSetJacobian()`

517:   Calling Sequence:
518: + snes   - the `SNES` context obtained from `SNESCreate()`
519: . u    - input vector
520: . Amat - (approximate) Jacobian matrix
521: . Pmat - matrix used to construct the preconditioner, often the same as `Amat`
522: - ctx  - [optional] user-defined context for matrix evaluation routine

524:   Level: beginner

526: .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESNGSFn`
527: S*/
528: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESJacobianFn(SNES snes, Vec u, Mat Amat, Mat Pmat, void *ctx);

530: /*S
531:   SNESNGSFn - A prototype of a `SNES` nonlinear Gauss-Seidel function that would be passed to `SNESSetNGS()`

533:   Calling Sequence:
534: + snes   - the `SNES` context obtained from `SNESCreate()`
535: . u    - the current solution, updated in place
536: . b    - the right-hand side vector (which may be `NULL`)
537: - ctx  - [optional] user-defined context for matrix evaluation routine

539:   Level: beginner

541: .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`
542: S*/
543: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESNGSFn(SNES snes, Vec u, Vec b, void *ctx);

545: /*S
546:   SNESUpdateFn - A prototype of a `SNES` update function that would be passed to `SNESSetUpdate()`

548:   Calling Sequence:
549: + snes - `SNES` context
550: - step - the current iteration index

552:   Level: advanced

554: .seealso: [](ch_snes), `SNES`, `SNESSetUpdate()`
555: S*/
556: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESUpdateFn(SNES snes, PetscInt step);

558: /* --------- Solving systems of nonlinear equations --------------- */
559: PETSC_EXTERN PetscErrorCode SNESSetFunction(SNES, Vec, SNESFunctionFn *, void *);
560: PETSC_EXTERN PetscErrorCode SNESGetFunction(SNES, Vec *, SNESFunctionFn **, void **);
561: PETSC_EXTERN PetscErrorCode SNESComputeFunction(SNES, Vec, Vec);
562: PETSC_EXTERN PetscErrorCode SNESComputeMFFunction(SNES, Vec, Vec);
563: PETSC_EXTERN PetscErrorCode SNESSetInitialFunction(SNES, Vec);

565: PETSC_EXTERN PetscErrorCode SNESSetJacobian(SNES, Mat, Mat, SNESJacobianFn *, void *);
566: PETSC_EXTERN PetscErrorCode SNESGetJacobian(SNES, Mat *, Mat *, SNESJacobianFn **, void **);
567: PETSC_EXTERN SNESFunctionFn SNESObjectiveComputeFunctionDefaultFD;
568: PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefault;
569: PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefaultColor;
570: PETSC_EXTERN PetscErrorCode SNESPruneJacobianColor(SNES, Mat, Mat);
571: PETSC_EXTERN PetscErrorCode SNESSetComputeInitialGuess(SNES, SNESInitialGuessFn *, void *);
572: PETSC_EXTERN PetscErrorCode SNESSetPicard(SNES, Vec, SNESFunctionFn *, Mat, Mat, SNESJacobianFn *, void *);
573: PETSC_EXTERN PetscErrorCode SNESGetPicard(SNES, Vec *, SNESFunctionFn **, Mat *, Mat *, SNESJacobianFn **, void **);
574: PETSC_EXTERN SNESFunctionFn SNESPicardComputeFunction;
575: PETSC_EXTERN SNESFunctionFn SNESPicardComputeMFFunction;
576: PETSC_EXTERN SNESJacobianFn SNESPicardComputeJacobian;

578: PETSC_EXTERN PetscErrorCode SNESSetObjective(SNES, SNESObjectiveFn *, void *);
579: PETSC_EXTERN PetscErrorCode SNESGetObjective(SNES, SNESObjectiveFn **, void **);
580: PETSC_EXTERN PetscErrorCode SNESComputeObjective(SNES, Vec, PetscReal *);

582: PETSC_EXTERN PetscErrorCode SNESSetUpdate(SNES, SNESUpdateFn *);

584: /*E
585:    SNESNormSchedule - Frequency with which the norm is computed during a nonliner solve

587:    Values:
588: +   `SNES_NORM_DEFAULT`            - use the default behavior for the current `SNESType`
589: .   `SNES_NORM_NONE`               - avoid all norm computations
590: .   `SNES_NORM_ALWAYS`             - compute the norms whenever possible
591: .   `SNES_NORM_INITIAL_ONLY`       - compute the norm only when the algorithm starts
592: .   `SNES_NORM_FINAL_ONLY`         - compute the norm only when the algorithm finishes
593: -   `SNES_NORM_INITIAL_FINAL_ONLY` - compute the norm at the start and end of the algorithm

595:    Level: advanced

597:    Notes:
598:    Support for these is highly dependent on the solver.

600:    Some options limit the convergence tests that can be used.

602:    The `SNES_NORM_NONE` option is most commonly used when the nonlinear solver is being used as a smoother, for example for `SNESFAS`

604:    This is primarily used to turn off extra norm and function computation
605:    when the solvers are composed.

607: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
608:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
609: E*/
610: typedef enum {
611:   SNES_NORM_DEFAULT            = -1,
612:   SNES_NORM_NONE               = 0,
613:   SNES_NORM_ALWAYS             = 1,
614:   SNES_NORM_INITIAL_ONLY       = 2,
615:   SNES_NORM_FINAL_ONLY         = 3,
616:   SNES_NORM_INITIAL_FINAL_ONLY = 4
617: } SNESNormSchedule;
618: PETSC_EXTERN const char *const *const SNESNormSchedules;

620: /*MC
621:    SNES_NORM_NONE - Don't compute function and its L2 norm when possible

623:    Level: advanced

625:    Note:
626:    This is most useful for stationary solvers with a fixed number of iterations used as smoothers.

628: .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_DEFAULT`
629: M*/

631: /*MC
632:    SNES_NORM_ALWAYS - Compute the function and its L2 norm at each iteration.

634:    Level: advanced

636:    Note:
637:    Most solvers will use this no matter what norm type is passed to them.

639: .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_NONE`
640: M*/

642: /*MC
643:    SNES_NORM_INITIAL_ONLY - Compute the function and its L2 at iteration 0, but do not update it.

645:    Level: advanced

647:    Notes:
648:    This method is useful in composed methods, when a true solution might actually be found before `SNESSolve()` is called.
649:    This option enables the solve to abort on the zeroth iteration if this is the case.

651:    For solvers that require the computation of the L2 norm of the function as part of the method, this merely cancels
652:    the norm computation at the last iteration (if possible).

654: .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_FINAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
655: M*/

657: /*MC
658:    SNES_NORM_FINAL_ONLY - Compute the function and its L2 norm on only the final iteration.

660:    Level: advanced

662:    Note:
663:    For solvers that require the computation of the L2 norm of the function as part of the method, behaves
664:    exactly as `SNES_NORM_DEFAULT`.  This method is useful when the function is gotten after `SNESSolve()` and
665:    used in subsequent computation for methods that do not need the norm computed during the rest of the
666:    solution procedure.

668: .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_INITIAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
669: M*/

671: /*MC
672:    SNES_NORM_INITIAL_FINAL_ONLY - Compute the function and its L2 norm on only the initial and final iterations.

674:    Level: advanced

676:    Note:
677:    This method combines the benefits of `SNES_NORM_INITIAL_ONLY` and `SNES_NORM_FINAL_ONLY`.

679: .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_SNES_NORM_INITIAL_ONLY`, `SNES_NORM_FINAL_ONLY`
680: M*/

682: PETSC_EXTERN PetscErrorCode SNESSetNormSchedule(SNES, SNESNormSchedule);
683: PETSC_EXTERN PetscErrorCode SNESGetNormSchedule(SNES, SNESNormSchedule *);
684: PETSC_EXTERN PetscErrorCode SNESSetFunctionNorm(SNES, PetscReal);
685: PETSC_EXTERN PetscErrorCode SNESGetFunctionNorm(SNES, PetscReal *);
686: PETSC_EXTERN PetscErrorCode SNESGetUpdateNorm(SNES, PetscReal *);
687: PETSC_EXTERN PetscErrorCode SNESGetSolutionNorm(SNES, PetscReal *);

689: /*E
690:    SNESFunctionType - Type of function computed

692:    Values:
693: +  `SNES_FUNCTION_DEFAULT`          - the default behavior for the current `SNESType`
694: .  `SNES_FUNCTION_UNPRECONDITIONED` - the original function provided
695: -  `SNES_FUNCTION_PRECONDITIONED`   - the modification of the function by the preconditioner

697:    Level: advanced

699:    Note:
700:    Support for these is dependent on the solver.

702: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
703:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
704: E*/
705: typedef enum {
706:   SNES_FUNCTION_DEFAULT          = -1,
707:   SNES_FUNCTION_UNPRECONDITIONED = 0,
708:   SNES_FUNCTION_PRECONDITIONED   = 1
709: } SNESFunctionType;
710: PETSC_EXTERN const char *const *const SNESFunctionTypes;

712: PETSC_EXTERN PetscErrorCode SNESSetFunctionType(SNES, SNESFunctionType);
713: PETSC_EXTERN PetscErrorCode SNESGetFunctionType(SNES, SNESFunctionType *);

715: PETSC_EXTERN PetscErrorCode SNESSetNGS(SNES, SNESNGSFn *, void *);
716: PETSC_EXTERN PetscErrorCode SNESGetNGS(SNES, SNESNGSFn **, void **);
717: PETSC_EXTERN PetscErrorCode SNESComputeNGS(SNES, Vec, Vec);

719: PETSC_EXTERN PetscErrorCode SNESNGSSetSweeps(SNES, PetscInt);
720: PETSC_EXTERN PetscErrorCode SNESNGSGetSweeps(SNES, PetscInt *);
721: PETSC_EXTERN PetscErrorCode SNESNGSSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt);
722: PETSC_EXTERN PetscErrorCode SNESNGSGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *);

724: PETSC_EXTERN PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES, PetscBool);
725: PETSC_EXTERN PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES, PetscBool *);

727: PETSC_EXTERN PetscErrorCode SNESShellGetContext(SNES, void *);
728: PETSC_EXTERN PetscErrorCode SNESShellSetContext(SNES, void *);
729: PETSC_EXTERN PetscErrorCode SNESShellSetSolve(SNES, PetscErrorCode (*)(SNES, Vec));

731: /* --------- Routines specifically for line search methods --------------- */

733: /*S
734:    SNESLineSearch - Abstract PETSc object that manages line-search operations for nonlinear solvers

736:    Level: beginner

738: .seealso: [](ch_snes), `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNES`
739: S*/
740: typedef struct _p_LineSearch *SNESLineSearch;

742: /*J
743:    SNESLineSearchType - String with the name of a PETSc line search method `SNESLineSearch`. Provides all the linesearches for the nonlinear solvers, `SNES`,
744:                         in PETSc.

746:    Values:
747: +  `SNESLINESEARCHBASIC`     - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
748: .  `SNESLINESEARCHBT`        - Backtracking line search over the L2 norm of the function or an objective function
749: .  `SNESLINESEARCHSECANT`    - Secant line search over the L2 norm of the function or an objective function
750: .  `SNESLINESEARCHCP`        - Critical point secant line search assuming $F(x) = \nabla G(x)$ for some unknown $G(x)$
751: .  `SNESLINESEARCHNLEQERR`   - Affine-covariant error-oriented linesearch
752: -  `SNESLINESEARCHBISECTION` - bisection line search for a root in the directional derivative
753: -  `SNESLINESEARCHSHELL`     - User provided `SNESLineSearch` implementation

755:    Level: beginner

757:    Note:
758:    Use `SNESLineSearchSetType()` or the options database key `-snes_linesearch_type` to set
759:    the specific line search algorithm to use with a given `SNES` object. Not all `SNESType` can utilize a line search.

761: .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetType()`, `SNES`
762: J*/
763: typedef const char *SNESLineSearchType;
764: #define SNESLINESEARCHBT        "bt"
765: #define SNESLINESEARCHNLEQERR   "nleqerr"
766: #define SNESLINESEARCHBASIC     "basic"
767: #define SNESLINESEARCHNONE      "none"
768: #define SNESLINESEARCHSECANT    "secant"
769: #define SNESLINESEARCHL2        PETSC_DEPRECATED_MACRO(3, 24, 0, "SNESLINESEARCHSECANT", ) "secant"
770: #define SNESLINESEARCHCP        "cp"
771: #define SNESLINESEARCHSHELL     "shell"
772: #define SNESLINESEARCHNCGLINEAR "ncglinear"
773: #define SNESLINESEARCHBISECTION "bisection"

775: PETSC_EXTERN PetscFunctionList SNESList;
776: PETSC_EXTERN PetscClassId      SNESLINESEARCH_CLASSID;
777: PETSC_EXTERN PetscFunctionList SNESLineSearchList;

779: #define SNES_LINESEARCH_ORDER_LINEAR    1
780: #define SNES_LINESEARCH_ORDER_QUADRATIC 2
781: #define SNES_LINESEARCH_ORDER_CUBIC     3

783: /*S
784:   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that projects a vector onto the VI bounds, passed to `SNESLineSearchSetVIFunctions()`

786:   Calling Sequence:
787: + snes  - `SNES` context
788: - u     - the vector to project to the bounds

790:   Level: advanced

792:   Note:
793:   The deprecated `SNESLineSearchVIProjectFunc` still works as a replacement for `SNESLineSearchVIProjectFn` *.

795: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
796: S*/
797: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode             SNESLineSearchVIProjectFn(SNES snes, Vec u);
798: PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVIProjectFn *SNESLineSearchVIProjectFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVIProjectFn*", );

800: /*S
801:   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that computes the norm of the active set variables in a vector in a VI solve,
802:   passed to `SNESLineSearchSetVIFunctions()`

804:   Calling Sequence:
805: + snes  - `SNES` context
806: . f     - the vector to compute the norm of
807: . u     - the current solution, entries that are on the VI bounds are ignored
808: - fnorm - the resulting norm

810:   Level: advanced

812:   Note:
813:   The deprecated `SNESLineSearchVINormFunc` still works as a replacement for `SNESLineSearchVINormFn` *.

815: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
816: S*/
817: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode          SNESLineSearchVINormFn(SNES snes, Vec f, Vec u, PetscReal *fnorm);
818: PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVINormFn *SNESLineSearchVINormFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVINormFnn*", );

820: /*S
821:   SNESLineSearchVIDirDerivFn - A prototype of a `SNES` function that computes the directional derivative considering the VI bounds, passed to `SNESLineSearchSetVIFunctions()`

823:   Calling Sequence:
824: + snes  - `SNES` context
825: . f     - the function vector to compute the directional derivative with
826: . u     - the current solution, entries that are on the VI bounds are ignored
827: . y     - the direction to compute the directional derivative
828: - fty   - the resulting directional derivative

830:   Level: advanced

832: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchVIProjectFn`, `SNESLineSearchVIProjectFn`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetVIFunctions()`
833: S*/
834: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchVIDirDerivFn(SNES snes, Vec f, Vec u, Vec y, PetscScalar *fty);

836: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchApplyFn(SNESLineSearch);
837: PETSC_EXTERN_TYPEDEF typedef SNESLineSearchApplyFn      *SNESLineSearchApplyFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
838: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchShellApplyFn(SNESLineSearch, void *);
839: PETSC_EXTERN_TYPEDEF typedef SNESLineSearchShellApplyFn *SNESLineSearchUserFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );

841: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate(MPI_Comm, SNESLineSearch *);
842: PETSC_EXTERN PetscErrorCode SNESLineSearchReset(SNESLineSearch);
843: PETSC_EXTERN PetscErrorCode SNESLineSearchView(SNESLineSearch, PetscViewer);
844: PETSC_EXTERN PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *);
845: PETSC_EXTERN PetscErrorCode SNESLineSearchGetType(SNESLineSearch, SNESLineSearchType *);
846: PETSC_EXTERN PetscErrorCode SNESLineSearchSetType(SNESLineSearch, SNESLineSearchType);
847: PETSC_EXTERN PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch);
848: PETSC_EXTERN PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch, PetscErrorCode (*)(SNES, Vec, Vec));
849: PETSC_EXTERN PetscErrorCode SNESLineSearchSetUp(SNESLineSearch);
850: PETSC_EXTERN PetscErrorCode SNESLineSearchApply(SNESLineSearch, Vec, Vec, PetscReal *, Vec);
851: PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch, Vec, Vec, PetscBool *);
852: PETSC_EXTERN PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *);
853: PETSC_EXTERN PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch, PetscInt);

855: /* set the functions for precheck and postcheck */

857: PETSC_EXTERN PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx);
858: PETSC_EXTERN PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);

860: PETSC_EXTERN PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx);
861: PETSC_EXTERN PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);

863: /* set the functions for VI-specific line search operations */

865: PETSC_EXTERN PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn *, SNESLineSearchVINormFn *, SNESLineSearchVIDirDerivFn *);
866: PETSC_EXTERN PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn **, SNESLineSearchVINormFn **, SNESLineSearchVIDirDerivFn **);

868: /* pointers to the associated SNES in order to be able to get the function evaluation out */
869: PETSC_EXTERN PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch, SNES);
870: PETSC_EXTERN PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch, SNES *);

872: /* set and get the parameters and vectors */
873: PETSC_EXTERN PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
874: PETSC_EXTERN PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscInt);

876: PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch, Vec, Vec, PetscBool *, void *);

878: PETSC_EXTERN PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch, PetscReal *);
879: PETSC_EXTERN PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch, PetscReal);

881: PETSC_EXTERN PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch, PetscReal *);
882: PETSC_EXTERN PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch, PetscReal);

884: PETSC_EXTERN PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch, PetscInt *);
885: PETSC_EXTERN PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch, PetscInt);

887: /*E
888:     SNESLineSearchReason - indication if the line search has succeeded or failed and why

890:   Values:
891: +  `SNES_LINESEARCH_SUCCEEDED`              - the line search succeeded
892: .  `SNES_LINESEARCH_FAILED_NANORINF`        - a not a number of infinity appeared in the computions
893: .  `SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN` - the function was evaluated outside of its domain, see `SNESSetFunctionDomainError()`
894: .  `SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN`- the objective function was evaluated outside of its domain, see `SNESSetObjectiveDomainError()`
895: .  `SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN` - the Jacobian was evaluated outside of its domain, see `SNESSetJacobianDomainError()`
896: .  `SNES_LINESEARCH_FAILED_REDUCT`          - the linear search failed to get the requested decrease in its norm or objective
897: .  `SNES_LINESEARCH_FAILED_USER`            - used by `SNESLINESEARCHNLEQERR` to indicate the user changed the search direction inappropriately
898: -  `SNES_LINESEARCH_FAILED_FUNCTION`        - indicates the maximum number of function evaluations allowed has been surpassed, `SNESConvergedReason` is also
899:                                               set to `SNES_DIVERGED_FUNCTION_COUNT`

901:    Level: intermediate

903:    Developer Note:
904:    Some of these reasons overlap with values of `SNESConvergedReason`. It is possibly a better design to have `SNESConvergedReaon` alone used also for indicating line
905:    search failures.

907: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`,
908:           `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
909: E*/
910: typedef enum {
911:   SNES_LINESEARCH_SUCCEEDED,
912:   SNES_LINESEARCH_FAILED_NANORINF,
913:   SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN,
914:   SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN,
915:   SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN,
916:   SNES_LINESEARCH_FAILED_REDUCT, /* INSUFFICIENT REDUCTION */
917:   SNES_LINESEARCH_FAILED_USER,
918:   SNES_LINESEARCH_FAILED_FUNCTION
919: } SNESLineSearchReason;

921: PETSC_EXTERN PetscErrorCode SNESLineSearchGetReason(SNESLineSearch, SNESLineSearchReason *);
922: PETSC_EXTERN PetscErrorCode SNESLineSearchSetReason(SNESLineSearch, SNESLineSearchReason);

924: PETSC_EXTERN PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch, Vec *, Vec *, Vec *, Vec *, Vec *);
925: PETSC_EXTERN PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch, Vec, Vec, Vec, Vec, Vec);

927: PETSC_EXTERN PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *);
928: PETSC_EXTERN PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch, PetscReal, PetscReal, PetscReal);
929: PETSC_EXTERN PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch);
930: PETSC_EXTERN PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch, PetscBool);

932: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitor(SNESLineSearch);
933: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, void *), void *, PetscCtxDestroyFn *);
934: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch, const char[], const char[], const char[], PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *));
935: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch);
936: PETSC_EXTERN PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch, PetscViewer);
937: PETSC_EXTERN PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch, PetscViewer *);
938: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch, PetscViewerAndFormat *);

940: PETSC_EXTERN PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch, const char[]);
941: PETSC_EXTERN PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch, const char *[]);

943: /* Shell interface functions */
944: PETSC_EXTERN PetscErrorCode SNESLineSearchShellSetApply(SNESLineSearch, SNESLineSearchShellApplyFn *, void *);
945: PETSC_EXTERN PetscErrorCode SNESLineSearchShellGetApply(SNESLineSearch, SNESLineSearchShellApplyFn **, void **);

947: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellSetApply()", ) static inline PetscErrorCode SNESLineSearchShellSetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn *f, void *ctx)
948: {
949:   return SNESLineSearchShellSetApply(ls, f, ctx);
950: }

952: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellGetApply()", ) static inline PetscErrorCode SNESLineSearchShellGetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn **f, void **ctx)
953: {
954:   return SNESLineSearchShellGetApply(ls, f, ctx);
955: }

957: /* BT interface functions */
958: PETSC_EXTERN PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch, PetscReal);
959: PETSC_EXTERN PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch, PetscReal *);

961: /*register line search types */
962: PETSC_EXTERN PetscErrorCode SNESLineSearchRegister(const char[], PetscErrorCode (*)(SNESLineSearch));

964: /* Routines for VI solver */
965: PETSC_EXTERN PetscErrorCode SNESVISetVariableBounds(SNES, Vec, Vec);
966: PETSC_EXTERN PetscErrorCode SNESVIGetVariableBounds(SNES, Vec *, Vec *);
967: PETSC_EXTERN PetscErrorCode SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
968: PETSC_EXTERN PetscErrorCode SNESVIGetInactiveSet(SNES, IS *);
969: PETSC_EXTERN PetscErrorCode SNESVIGetActiveSetIS(SNES, Vec, Vec, IS *);
970: PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES, Vec, Vec, PetscReal *);
971: PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFtY(SNES, Vec, Vec, Vec, PetscScalar *);
972: PETSC_EXTERN PetscErrorCode SNESVISetRedundancyCheck(SNES, PetscErrorCode (*)(SNES, IS, IS *, void *), void *);
973: PETSC_EXTERN PetscErrorCode SNESVIComputeMeritFunction(Vec, PetscReal *, PetscReal *);
974: PETSC_EXTERN PetscErrorCode SNESVIComputeFunction(SNES, Vec, Vec, void *);
975: PETSC_EXTERN PetscErrorCode DMSetVI(DM, IS);
976: PETSC_EXTERN PetscErrorCode DMDestroyVI(DM);

978: PETSC_EXTERN PetscErrorCode SNESTestLocalMin(SNES);

980: /* Should this routine be private? */
981: PETSC_EXTERN PetscErrorCode SNESComputeJacobian(SNES, Vec, Mat, Mat);
982: PETSC_EXTERN PetscErrorCode SNESTestJacobian(SNES, PetscReal *, PetscReal *);
983: PETSC_EXTERN PetscErrorCode SNESTestFunction(SNES);

985: PETSC_EXTERN PetscErrorCode SNESSetDM(SNES, DM);
986: PETSC_EXTERN PetscErrorCode SNESGetDM(SNES, DM *);
987: PETSC_EXTERN PetscErrorCode SNESSetNPC(SNES, SNES);
988: PETSC_EXTERN PetscErrorCode SNESGetNPC(SNES, SNES *);
989: PETSC_EXTERN PetscErrorCode SNESHasNPC(SNES, PetscBool *);
990: PETSC_EXTERN PetscErrorCode SNESApplyNPC(SNES, Vec, Vec, Vec);
991: PETSC_EXTERN PetscErrorCode SNESGetNPCFunction(SNES, Vec, PetscReal *);
992: PETSC_EXTERN PetscErrorCode SNESComputeFunctionDefaultNPC(SNES, Vec, Vec);
993: PETSC_EXTERN PetscErrorCode SNESSetNPCSide(SNES, PCSide);
994: PETSC_EXTERN PetscErrorCode SNESGetNPCSide(SNES, PCSide *);
995: PETSC_EXTERN PetscErrorCode SNESSetLineSearch(SNES, SNESLineSearch);
996: PETSC_EXTERN PetscErrorCode SNESGetLineSearch(SNES, SNESLineSearch *);

998: PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESGetLineSearch()", ) static inline PetscErrorCode SNESGetSNESLineSearch(SNES snes, SNESLineSearch *ls)
999: {
1000:   return SNESGetLineSearch(snes, ls);
1001: }
1002: PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESSetLineSearch()", ) static inline PetscErrorCode SNESSetSNESLineSearch(SNES snes, SNESLineSearch ls)
1003: {
1004:   return SNESSetLineSearch(snes, ls);
1005: }

1007: PETSC_EXTERN PetscErrorCode SNESSetUpMatrices(SNES);
1008: PETSC_EXTERN PetscErrorCode DMSNESSetFunction(DM, SNESFunctionFn *, void *);
1009: PETSC_EXTERN PetscErrorCode DMSNESGetFunction(DM, SNESFunctionFn **, void **);
1010: PETSC_EXTERN PetscErrorCode DMSNESSetFunctionContextDestroy(DM, PetscCtxDestroyFn *);
1011: PETSC_EXTERN PetscErrorCode DMSNESSetMFFunction(DM, SNESFunctionFn *, void *);
1012: PETSC_EXTERN PetscErrorCode DMSNESSetNGS(DM, SNESNGSFn *, void *);
1013: PETSC_EXTERN PetscErrorCode DMSNESGetNGS(DM, SNESNGSFn **, void **);
1014: PETSC_EXTERN PetscErrorCode DMSNESSetJacobian(DM, SNESJacobianFn *, void *);
1015: PETSC_EXTERN PetscErrorCode DMSNESGetJacobian(DM, SNESJacobianFn **, void **);
1016: PETSC_EXTERN PetscErrorCode DMSNESSetJacobianContextDestroy(DM, PetscCtxDestroyFn *);
1017: PETSC_EXTERN PetscErrorCode DMSNESSetPicard(DM, SNESFunctionFn *, SNESJacobianFn *, void *);
1018: PETSC_EXTERN PetscErrorCode DMSNESGetPicard(DM, SNESFunctionFn **, SNESJacobianFn **, void **);
1019: PETSC_EXTERN PetscErrorCode DMSNESSetObjective(DM, SNESObjectiveFn *, void *);
1020: PETSC_EXTERN PetscErrorCode DMSNESGetObjective(DM, SNESObjectiveFn **, void **);
1021: PETSC_EXTERN PetscErrorCode DMCopyDMSNES(DM, DM);

1023: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionFn(DMDALocalInfo *, void *, void *, void *);
1024: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianFn(DMDALocalInfo *, void *, Mat, Mat, void *);
1025: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveFn(DMDALocalInfo *, void *, PetscReal *, void *);

1027: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionVecFn(DMDALocalInfo *, Vec, Vec, void *);
1028: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianVecFn(DMDALocalInfo *, Vec, Mat, Mat, void *);
1029: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveVecFn(DMDALocalInfo *, Vec, PetscReal *, void *);

1031: PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocal(DM, InsertMode, DMDASNESFunctionFn *, void *);
1032: PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocal(DM, DMDASNESJacobianFn *, void *);
1033: PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocal(DM, DMDASNESObjectiveFn *, void *);
1034: PETSC_EXTERN PetscErrorCode DMDASNESSetPicardLocal(DM, InsertMode, DMDASNESFunctionFn *, DMDASNESJacobianFn, void *);

1036: PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocalVec(DM, InsertMode, DMDASNESFunctionVecFn *, void *);
1037: PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocalVec(DM, DMDASNESJacobianVecFn *, void *);
1038: PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocalVec(DM, DMDASNESObjectiveVecFn *, void *);

1040: PETSC_EXTERN PetscErrorCode DMSNESSetBoundaryLocal(DM, PetscErrorCode (*)(DM, Vec, void *), void *);
1041: PETSC_EXTERN PetscErrorCode DMSNESSetObjectiveLocal(DM, PetscErrorCode (*)(DM, Vec, PetscReal *, void *), void *);
1042: PETSC_EXTERN PetscErrorCode DMSNESSetFunctionLocal(DM, PetscErrorCode (*)(DM, Vec, Vec, void *), void *);
1043: PETSC_EXTERN PetscErrorCode DMSNESSetJacobianLocal(DM, PetscErrorCode (*)(DM, Vec, Mat, Mat, void *), void *);
1044: PETSC_EXTERN PetscErrorCode DMSNESGetBoundaryLocal(DM, PetscErrorCode (**)(DM, Vec, void *), void **);
1045: PETSC_EXTERN PetscErrorCode DMSNESGetObjectiveLocal(DM, PetscErrorCode (**)(DM, Vec, PetscReal *, void *), void **);
1046: PETSC_EXTERN PetscErrorCode DMSNESGetFunctionLocal(DM, PetscErrorCode (**)(DM, Vec, Vec, void *), void **);
1047: PETSC_EXTERN PetscErrorCode DMSNESGetJacobianLocal(DM, PetscErrorCode (**)(DM, Vec, Mat, Mat, void *), void **);

1049: /* Routines for Multiblock solver */
1050: PETSC_EXTERN PetscErrorCode SNESMultiblockSetFields(SNES, const char[], PetscInt, const PetscInt *);
1051: PETSC_EXTERN PetscErrorCode SNESMultiblockSetIS(SNES, const char[], IS);
1052: PETSC_EXTERN PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
1053: PETSC_EXTERN PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);
1054: PETSC_EXTERN PetscErrorCode SNESMultiblockGetSubSNES(SNES, PetscInt *, SNES *[]);

1056: /*J
1057:    SNESMSType - String with the name of a PETSc `SNESMS` method.

1059:    Level: intermediate

1061: .seealso: [](ch_snes), `SNESMS`, `SNESMSGetType()`, `SNESMSSetType()`, `SNES`
1062: J*/
1063: typedef const char *SNESMSType;
1064: #define SNESMSM62       "m62"
1065: #define SNESMSEULER     "euler"
1066: #define SNESMSJAMESON83 "jameson83"
1067: #define SNESMSVLTP11    "vltp11"
1068: #define SNESMSVLTP21    "vltp21"
1069: #define SNESMSVLTP31    "vltp31"
1070: #define SNESMSVLTP41    "vltp41"
1071: #define SNESMSVLTP51    "vltp51"
1072: #define SNESMSVLTP61    "vltp61"

1074: PETSC_EXTERN PetscErrorCode SNESMSRegister(SNESMSType, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[]);
1075: PETSC_EXTERN PetscErrorCode SNESMSRegisterAll(void);
1076: PETSC_EXTERN PetscErrorCode SNESMSGetType(SNES, SNESMSType *);
1077: PETSC_EXTERN PetscErrorCode SNESMSSetType(SNES, SNESMSType);
1078: PETSC_EXTERN PetscErrorCode SNESMSGetDamping(SNES, PetscReal *);
1079: PETSC_EXTERN PetscErrorCode SNESMSSetDamping(SNES, PetscReal);
1080: PETSC_EXTERN PetscErrorCode SNESMSFinalizePackage(void);
1081: PETSC_EXTERN PetscErrorCode SNESMSInitializePackage(void);
1082: PETSC_EXTERN PetscErrorCode SNESMSRegisterDestroy(void);

1084: /*MC
1085:    SNESNGMRESRestartType - the restart approach used by `SNESNGMRES`

1087:   Values:
1088: +   `SNES_NGMRES_RESTART_NONE`       - never restart
1089: .   `SNES_NGMRES_RESTART_DIFFERENCE` - restart based upon difference criteria
1090: -   `SNES_NGMRES_RESTART_PERIODIC`   - restart after a fixed number of iterations

1092:   Options Database Keys:
1093: + -snes_ngmres_restart_type <difference,periodic,none> - set the restart type
1094: - -snes_ngmres_restart <30>                            - sets the number of iterations before restart for periodic

1096:    Level: intermediate

1098: .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1099:           `SNESNGMRESGetRestartType()`, `SNESNGMRESSelectType`
1100: M*/
1101: typedef enum {
1102:   SNES_NGMRES_RESTART_NONE       = 0,
1103:   SNES_NGMRES_RESTART_PERIODIC   = 1,
1104:   SNES_NGMRES_RESTART_DIFFERENCE = 2
1105: } SNESNGMRESRestartType;
1106: PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];

1108: /*MC
1109:    SNESNGMRESSelectType - the approach used by `SNESNGMRES` to determine how the candidate solution and
1110:   combined solution are used to create the next iterate.

1112:    Values:
1113: +   `SNES_NGMRES_SELECT_NONE`       - choose the combined solution all the time
1114: .   `SNES_NGMRES_SELECT_DIFFERENCE` - choose based upon the selection criteria
1115: -   `SNES_NGMRES_SELECT_LINESEARCH` - choose based upon line search combination

1117:   Options Database Key:
1118: . -snes_ngmres_select_type<difference,none,linesearch> - select type

1120:    Level: intermediate

1122: .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1123:           `SNESNGMRESGetRestartType()`, `SNESNGMRESRestartType`
1124: M*/
1125: typedef enum {
1126:   SNES_NGMRES_SELECT_NONE       = 0,
1127:   SNES_NGMRES_SELECT_DIFFERENCE = 1,
1128:   SNES_NGMRES_SELECT_LINESEARCH = 2
1129: } SNESNGMRESSelectType;
1130: PETSC_EXTERN const char *const SNESNGMRESSelectTypes[];

1132: PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartType(SNES, SNESNGMRESRestartType);
1133: PETSC_EXTERN PetscErrorCode SNESNGMRESSetSelectType(SNES, SNESNGMRESSelectType);
1134: PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartFmRise(SNES, PetscBool);
1135: PETSC_EXTERN PetscErrorCode SNESNGMRESGetRestartFmRise(SNES, PetscBool *);

1137: /*MC
1138:    SNESNCGType - the conjugate update approach for `SNESNCG`

1140:    Values:
1141: +   `SNES_NCG_FR`  - Fletcher-Reeves update
1142: .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update, the default and the only one that tolerates generalized search directions
1143: .   `SNES_NCG_HS`  - Hestenes-Steifel update
1144: .   `SNES_NCG_DY`  - Dai-Yuan update
1145: -   `SNES_NCG_CD`  - Conjugate Descent update

1147:   Options Database Key:
1148: . -snes_ncg_type<fr,prp,hs,dy,cd> - select type

1150:    Level: intermediate

1152: .seealso: `SNES, `SNESNCG`, `SNESNCGSetType()`
1153: M*/
1154: typedef enum {
1155:   SNES_NCG_FR  = 0,
1156:   SNES_NCG_PRP = 1,
1157:   SNES_NCG_HS  = 2,
1158:   SNES_NCG_DY  = 3,
1159:   SNES_NCG_CD  = 4
1160: } SNESNCGType;
1161: PETSC_EXTERN const char *const SNESNCGTypes[];

1163: PETSC_EXTERN PetscErrorCode SNESNCGSetType(SNES, SNESNCGType);

1165: /*MC
1166:    SNESQNScaleType - the scaling type used by `SNESQN`

1168:    Values:
1169: +   `SNES_QN_SCALE_NONE`     - don't scale the problem
1170: .   `SNES_QN_SCALE_SCALAR`   - use Shanno scaling
1171: .   `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
1172: -   `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with `SNESSetJacobian()`
1173:                                computed at the first iteration of `SNESQN` and at ever restart.

1175:     Options Database Key:
1176: . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type

1178:    Level: intermediate

1180: .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNRestartType`
1181: M*/
1182: typedef enum {
1183:   SNES_QN_SCALE_DEFAULT  = 0,
1184:   SNES_QN_SCALE_NONE     = 1,
1185:   SNES_QN_SCALE_SCALAR   = 2,
1186:   SNES_QN_SCALE_DIAGONAL = 3,
1187:   SNES_QN_SCALE_JACOBIAN = 4
1188: } SNESQNScaleType;
1189: PETSC_EXTERN const char *const SNESQNScaleTypes[];

1191: /*MC
1192:    SNESQNRestartType - the restart approached used by `SNESQN`

1194:    Values:
1195: +   `SNES_QN_RESTART_NONE`     - never restart
1196: .   `SNES_QN_RESTART_POWELL`   - restart based upon descent criteria
1197: -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations

1199:   Options Database Keys:
1200: + -snes_qn_restart_type <powell,periodic,none> - set the restart type
1201: - -snes_qn_m <m>                               - sets the number of stored updates and the restart period for periodic

1203:    Level: intermediate

1205: .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNScaleType`
1206: M*/
1207: typedef enum {
1208:   SNES_QN_RESTART_DEFAULT  = 0,
1209:   SNES_QN_RESTART_NONE     = 1,
1210:   SNES_QN_RESTART_POWELL   = 2,
1211:   SNES_QN_RESTART_PERIODIC = 3
1212: } SNESQNRestartType;
1213: PETSC_EXTERN const char *const SNESQNRestartTypes[];

1215: /*MC
1216:    SNESQNType - the type used by `SNESQN`

1218:   Values:
1219: +   `SNES_QN_LBFGS`      - LBFGS variant
1220: .   `SNES_QN_BROYDEN`    - Broyden variant
1221: -   `SNES_QN_BADBROYDEN` - Bad Broyden variant

1223:   Options Database Key:
1224: . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type

1226:    Level: intermediate

1228: .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNSetType()`, `SNESQNScaleType`, `SNESQNRestartType`, `SNESQNSetRestartType()`
1229: M*/
1230: typedef enum {
1231:   SNES_QN_LBFGS      = 0,
1232:   SNES_QN_BROYDEN    = 1,
1233:   SNES_QN_BADBROYDEN = 2
1234: } SNESQNType;
1235: PETSC_EXTERN const char *const SNESQNTypes[];

1237: PETSC_EXTERN PetscErrorCode SNESQNSetType(SNES, SNESQNType);
1238: PETSC_EXTERN PetscErrorCode SNESQNSetScaleType(SNES, SNESQNScaleType);
1239: PETSC_EXTERN PetscErrorCode SNESQNSetRestartType(SNES, SNESQNRestartType);

1241: PETSC_EXTERN PetscErrorCode SNESNASMGetType(SNES, PCASMType *);
1242: PETSC_EXTERN PetscErrorCode SNESNASMSetType(SNES, PCASMType);
1243: PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomains(SNES, PetscInt *, SNES *[], VecScatter *[], VecScatter *[], VecScatter *[]);
1244: PETSC_EXTERN PetscErrorCode SNESNASMSetSubdomains(SNES, PetscInt, SNES[], VecScatter[], VecScatter[], VecScatter[]);
1245: PETSC_EXTERN PetscErrorCode SNESNASMSetDamping(SNES, PetscReal);
1246: PETSC_EXTERN PetscErrorCode SNESNASMGetDamping(SNES, PetscReal *);
1247: PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomainVecs(SNES, PetscInt *, Vec *[], Vec *[], Vec *[], Vec *[]);
1248: PETSC_EXTERN PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES, PetscBool);
1249: PETSC_EXTERN PetscErrorCode SNESNASMGetSNES(SNES, PetscInt, SNES *);
1250: PETSC_EXTERN PetscErrorCode SNESNASMGetNumber(SNES, PetscInt *);
1251: PETSC_EXTERN PetscErrorCode SNESNASMSetWeight(SNES, Vec);

1253: /*E
1254:   SNESCompositeType - Determines how two or more preconditioners are composed with the `SNESType` of `SNESCOMPOSITE`

1256:   Values:
1257: + `SNES_COMPOSITE_ADDITIVE`        - results from application of all preconditioners are added together
1258: . `SNES_COMPOSITE_MULTIPLICATIVE`  - preconditioners are applied sequentially to the residual freshly
1259:                                      computed after the previous preconditioner application
1260: - `SNES_COMPOSITE_ADDITIVEOPTIMAL` - uses a linear combination of the solutions obtained with each preconditioner that approximately minimize the function
1261:                                      value at the new iteration.

1263:    Level: beginner

1265: .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `PCCompositeType`
1266: E*/
1267: typedef enum {
1268:   SNES_COMPOSITE_ADDITIVE,
1269:   SNES_COMPOSITE_MULTIPLICATIVE,
1270:   SNES_COMPOSITE_ADDITIVEOPTIMAL
1271: } SNESCompositeType;
1272: PETSC_EXTERN const char *const SNESCompositeTypes[];

1274: PETSC_EXTERN PetscErrorCode SNESCompositeSetType(SNES, SNESCompositeType);
1275: PETSC_EXTERN PetscErrorCode SNESCompositeAddSNES(SNES, SNESType);
1276: PETSC_EXTERN PetscErrorCode SNESCompositeGetSNES(SNES, PetscInt, SNES *);
1277: PETSC_EXTERN PetscErrorCode SNESCompositeGetNumber(SNES, PetscInt *);
1278: PETSC_EXTERN PetscErrorCode SNESCompositeSetDamping(SNES, PetscInt, PetscReal);

1280: PETSC_EXTERN PetscErrorCode SNESPatchSetDiscretisationInfo(SNES, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
1281: PETSC_EXTERN PetscErrorCode SNESPatchSetComputeOperator(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
1282: PETSC_EXTERN PetscErrorCode SNESPatchSetComputeFunction(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
1283: PETSC_EXTERN PetscErrorCode SNESPatchSetConstructType(SNES, PCPatchConstructType, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *);
1284: PETSC_EXTERN PetscErrorCode SNESPatchSetCellNumbering(SNES, PetscSection);

1286: /*E
1287:     SNESFASType - Determines the type of nonlinear multigrid method that is run.

1289:    Values:
1290: +  `SNES_FAS_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `SNESFASSetCycles()`
1291: .  `SNES_FAS_ADDITIVE`                 - additive FAS cycle
1292: .  `SNES_FAS_FULL`                     - full FAS cycle
1293: -  `SNES_FAS_KASKADE`                  - Kaskade FAS cycle

1295:    Level: beginner

1297: .seealso: [](ch_snes), `SNESFAS`, `PCMGSetType()`, `PCMGType`
1298: E*/
1299: typedef enum {
1300:   SNES_FAS_MULTIPLICATIVE,
1301:   SNES_FAS_ADDITIVE,
1302:   SNES_FAS_FULL,
1303:   SNES_FAS_KASKADE
1304: } SNESFASType;
1305: PETSC_EXTERN const char *const SNESFASTypes[];

1307: /* called on the finest level FAS instance*/
1308: PETSC_EXTERN PetscErrorCode SNESFASSetType(SNES, SNESFASType);
1309: PETSC_EXTERN PetscErrorCode SNESFASGetType(SNES, SNESFASType *);
1310: PETSC_EXTERN PetscErrorCode SNESFASSetLevels(SNES, PetscInt, MPI_Comm *);
1311: PETSC_EXTERN PetscErrorCode SNESFASGetLevels(SNES, PetscInt *);
1312: PETSC_EXTERN PetscErrorCode SNESFASGetCycleSNES(SNES, PetscInt, SNES *);
1313: PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothUp(SNES, PetscInt);
1314: PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothDown(SNES, PetscInt);
1315: PETSC_EXTERN PetscErrorCode SNESFASSetCycles(SNES, PetscInt);
1316: PETSC_EXTERN PetscErrorCode SNESFASSetMonitor(SNES, PetscViewerAndFormat *, PetscBool);
1317: PETSC_EXTERN PetscErrorCode SNESFASSetLog(SNES, PetscBool);

1319: PETSC_EXTERN PetscErrorCode SNESFASSetGalerkin(SNES, PetscBool);
1320: PETSC_EXTERN PetscErrorCode SNESFASGetGalerkin(SNES, PetscBool *);
1321: PETSC_EXTERN PetscErrorCode SNESFASGalerkinFunctionDefault(SNES, Vec, Vec, void *);

1323: /* called on any level -- "Cycle" FAS instance */
1324: PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmoother(SNES, SNES *);
1325: PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherUp(SNES, SNES *);
1326: PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherDown(SNES, SNES *);
1327: PETSC_EXTERN PetscErrorCode SNESFASCycleGetCorrection(SNES, SNES *);
1328: PETSC_EXTERN PetscErrorCode SNESFASCycleGetInterpolation(SNES, Mat *);
1329: PETSC_EXTERN PetscErrorCode SNESFASCycleGetRestriction(SNES, Mat *);
1330: PETSC_EXTERN PetscErrorCode SNESFASCycleGetInjection(SNES, Mat *);
1331: PETSC_EXTERN PetscErrorCode SNESFASCycleGetRScale(SNES, Vec *);
1332: PETSC_EXTERN PetscErrorCode SNESFASCycleSetCycles(SNES, PetscInt);
1333: PETSC_EXTERN PetscErrorCode SNESFASCycleIsFine(SNES, PetscBool *);

1335: /* called on the (outer) finest level FAS to set/get parameters on any level instance */
1336: PETSC_EXTERN PetscErrorCode SNESFASSetInterpolation(SNES, PetscInt, Mat);
1337: PETSC_EXTERN PetscErrorCode SNESFASGetInterpolation(SNES, PetscInt, Mat *);
1338: PETSC_EXTERN PetscErrorCode SNESFASSetRestriction(SNES, PetscInt, Mat);
1339: PETSC_EXTERN PetscErrorCode SNESFASGetRestriction(SNES, PetscInt, Mat *);
1340: PETSC_EXTERN PetscErrorCode SNESFASSetInjection(SNES, PetscInt, Mat);
1341: PETSC_EXTERN PetscErrorCode SNESFASGetInjection(SNES, PetscInt, Mat *);
1342: PETSC_EXTERN PetscErrorCode SNESFASSetRScale(SNES, PetscInt, Vec);
1343: PETSC_EXTERN PetscErrorCode SNESFASGetRScale(SNES, PetscInt, Vec *);
1344: PETSC_EXTERN PetscErrorCode SNESFASSetContinuation(SNES, PetscBool);

1346: PETSC_EXTERN PetscErrorCode SNESFASGetSmoother(SNES, PetscInt, SNES *);
1347: PETSC_EXTERN PetscErrorCode SNESFASGetSmootherUp(SNES, PetscInt, SNES *);
1348: PETSC_EXTERN PetscErrorCode SNESFASGetSmootherDown(SNES, PetscInt, SNES *);
1349: PETSC_EXTERN PetscErrorCode SNESFASGetCoarseSolve(SNES, SNES *);

1351: /* parameters for full FAS */
1352: PETSC_EXTERN PetscErrorCode SNESFASFullSetDownSweep(SNES, PetscBool);
1353: PETSC_EXTERN PetscErrorCode SNESFASCreateCoarseVec(SNES, Vec *);
1354: PETSC_EXTERN PetscErrorCode SNESFASRestrict(SNES, Vec, Vec);
1355: PETSC_EXTERN PetscErrorCode SNESFASFullSetTotal(SNES, PetscBool);
1356: PETSC_EXTERN PetscErrorCode SNESFASFullGetTotal(SNES, PetscBool *);

1358: PETSC_EXTERN PetscErrorCode DMPlexSetSNESVariableBounds(DM, SNES);
1359: PETSC_EXTERN PetscErrorCode DMSNESCheckDiscretization(SNES, DM, PetscReal, Vec, PetscReal, PetscReal[]);
1360: PETSC_EXTERN PetscErrorCode DMSNESCheckResidual(SNES, DM, Vec, PetscReal, PetscReal *);
1361: PETSC_EXTERN PetscErrorCode DMSNESCheckJacobian(SNES, DM, Vec, PetscReal, PetscBool *, PetscReal *);
1362: PETSC_EXTERN PetscErrorCode DMSNESCheckFromOptions(SNES, Vec);
1363: PETSC_EXTERN PetscErrorCode DMSNESComputeJacobianAction(DM, Vec, Vec, Vec, void *);
1364: PETSC_EXTERN PetscErrorCode DMSNESCreateJacobianMF(DM, Vec, void *, Mat *);

1366: PETSC_EXTERN PetscErrorCode SNESNewtonALSetFunction(SNES, SNESFunctionFn *, void *ctx);
1367: PETSC_EXTERN PetscErrorCode SNESNewtonALGetFunction(SNES, SNESFunctionFn **, void **ctx);
1368: PETSC_EXTERN PetscErrorCode SNESNewtonALComputeFunction(SNES, Vec, Vec);
1369: PETSC_EXTERN PetscErrorCode SNESNewtonALGetLoadParameter(SNES, PetscReal *);

1371: /*MC
1372:    SNESNewtonALCorrectionType - the approach used by `SNESNEWTONAL` to determine
1373:    the correction to the current increment. While the exact correction satisfies
1374:    the constraint surface at every iteration, it also requires solving a quadratic
1375:    equation which may not have real roots. Conversely, the normal correction is more
1376:    efficient and always yields a real correction and is the default.

1378:    Values:
1379: +   `SNES_NEWTONAL_CORRECTION_EXACT` - choose the correction which exactly satisfies the constraint
1380: -   `SNES_NEWTONAL_CORRECTION_NORMAL` - choose the correction in the updated normal hyper-surface to the constraint surface

1382:    Options Database Key:
1383: . -snes_newtonal_correction_type <exact> - select type from <exact,normal>

1385:    Level: intermediate

1387: .seealso: `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetCorrectionType()`
1388: M*/
1389: typedef enum {
1390:   SNES_NEWTONAL_CORRECTION_EXACT  = 0,
1391:   SNES_NEWTONAL_CORRECTION_NORMAL = 1,
1392: } SNESNewtonALCorrectionType;
1393: PETSC_EXTERN const char *const SNESNewtonALCorrectionTypes[];

1395: PETSC_EXTERN PetscErrorCode SNESNewtonALSetCorrectionType(SNES, SNESNewtonALCorrectionType);