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 SNESGetFunctionDomainError(SNES, PetscBool *);
249: PETSC_EXTERN PetscErrorCode SNESGetJacobianDomainError(SNES, PetscBool *);
250: PETSC_EXTERN PetscErrorCode SNESSetJacobianDomainError(SNES);
251: PETSC_EXTERN PetscErrorCode SNESSetCheckJacobianDomainError(SNES, PetscBool);
252: PETSC_EXTERN PetscErrorCode SNESGetCheckJacobianDomainError(SNES, PetscBool *);

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

258:    Values:
259: +  `SNES_CONVERGED_FNORM_ABS`      - 2-norm(F) <= abstol
260: .  `SNES_CONVERGED_FNORM_RELATIVE` - 2-norm(F) <= rtol*2-norm(F(x_0)) where x_0 is the initial guess
261: .  `SNES_CONVERGED_SNORM_RELATIVE` - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
262: .  `SNES_CONVERGED_USER`           - The user has indicated convergence for an arbitrary reason
263: .  `SNES_DIVERGED_FUNCTION_COUNT`  - The user provided function has been called more times than the maximum set in `SNESSetTolerances()`
264: .  `SNES_DIVERGED_DTOL`            - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
265: .  `SNES_DIVERGED_FNORM_NAN`       - the 2-norm of the current function evaluation is not-a-number (NaN), this
266:                                      is usually caused by a division of 0 by 0.
267: .  `SNES_DIVERGED_MAX_IT`          - `SNESSolve()` has reached the maximum number of iterations requested
268: .  `SNES_DIVERGED_LINE_SEARCH`     - The line search has failed. This only occurs for `SNES` solvers that use a line search
269: .  `SNES_DIVERGED_LOCAL_MIN`       - the algorithm seems to have stagnated at a local minimum that is not zero.
270: .  `SNES_DIVERGED_USER`            - The user has indicated divergence for an arbitrary reason
271: -  `SNES_CONVERGED_ITERATING       - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`

273:    Level: beginner

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

280:    `SNES_DIVERGED_LOCAL_MIN` can only occur when using a `SNES` solver that uses a line search (`SNESLineSearch`).
281:    The line search wants to minimize Q(alpha) = 1/2 || F(x + alpha s) ||^2_2  this occurs
282:    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
283:    you get Q'(alpha) = -F(x)^T F'(x)^(-1)^T F'(x+alpha s)F(x+alpha s); when alpha = 0
284:    Q'(0) = - ||F(x)||^2_2 which is always NEGATIVE if F'(x) is invertible. This means the Newton
285:    direction is a descent direction and the line search should succeed if alpha is small enough.

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

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

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

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

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

332:   SNES_CONVERGED_ITERATING = 0
333: } SNESConvergedReason;
334: PETSC_EXTERN const char *const *SNESConvergedReasons;

336: /*MC
337:    SNES_CONVERGED_FNORM_ABS - 2-norm(F) <= abstol

339:    Level: beginner

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

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

347:    Level: beginner

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

352: /*MC
353:   SNES_CONVERGED_SNORM_RELATIVE - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
354:   solution and stol is the 4th argument to `SNESSetTolerances()`

356:   Options Database Key:
357:   -snes_stol <stol> - the step tolerance

359:    Level: beginner

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

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

368:    Level: beginner

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

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

376:    Level: beginner

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

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

385:    Level: beginner

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

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

393:    Level: beginner

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

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

401:    Level: beginner

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

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

410:    Level: beginner

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

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

418:    Level: beginner

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

423: PETSC_EXTERN PetscErrorCode SNESSetConvergenceTest(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
424: PETSC_EXTERN PetscErrorCode SNESConvergedDefault(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
425: PETSC_EXTERN PetscErrorCode SNESConvergedSkip(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
426: PETSC_EXTERN PetscErrorCode SNESConvergedCorrectPressure(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
427: PETSC_EXTERN PetscErrorCode SNESGetConvergedReason(SNES, SNESConvergedReason *);
428: PETSC_EXTERN PetscErrorCode SNESGetConvergedReasonString(SNES, const char **);
429: PETSC_EXTERN PetscErrorCode SNESSetConvergedReason(SNES, SNESConvergedReason);

431: PETSC_DEPRECATED_FUNCTION(3, 5, 0, "SNESConvergedSkip()", ) static inline void SNESSkipConverged(void)
432: { /* never called */
433: }
434: #define SNESSkipConverged (SNESSkipConverged, SNESConvergedSkip)

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

439:   Calling Sequence:
440: + snes  - `SNES` context
441: . u   - output vector to contain initial guess
442: - ctx - [optional] user-defined function context

444:   Level: beginner

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

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

453:   Calling Sequence:
454: + snes  - `SNES` context
455: . u   - input vector
456: . F   - function vector
457: - ctx - [optional] user-defined function context

459:   Level: beginner

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

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

468:   Calling Sequence:
469: + snes  - `SNES` context
470: . u   - input vector
471: . o   - output value
472: - ctx - [optional] user-defined function context

474:   Level: beginner

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

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

483:   Calling Sequence:
484: + snes   - the `SNES` context obtained from `SNESCreate()`
485: . u    - input vector
486: . Amat - (approximate) Jacobian matrix
487: . Pmat - matrix used to construct the preconditioner, often the same as `Amat`
488: - ctx  - [optional] user-defined context for matrix evaluation routine

490:   Level: beginner

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

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

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

505:   Level: beginner

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

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

514:   Calling Sequence:
515: + snes - `SNES` context
516: - step - the current iteration index

518:   Level: advanced

520: .seealso: [](ch_snes), `SNES`, `SNESSetUpdate()`
521: S*/
522: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESUpdateFn(SNES snes, PetscInt step);

524: /* --------- Solving systems of nonlinear equations --------------- */
525: PETSC_EXTERN PetscErrorCode SNESSetFunction(SNES, Vec, SNESFunctionFn *, void *);
526: PETSC_EXTERN PetscErrorCode SNESGetFunction(SNES, Vec *, SNESFunctionFn **, void **);
527: PETSC_EXTERN PetscErrorCode SNESComputeFunction(SNES, Vec, Vec);
528: PETSC_EXTERN PetscErrorCode SNESComputeMFFunction(SNES, Vec, Vec);
529: PETSC_EXTERN PetscErrorCode SNESSetInitialFunction(SNES, Vec);

531: PETSC_EXTERN PetscErrorCode SNESSetJacobian(SNES, Mat, Mat, SNESJacobianFn *, void *);
532: PETSC_EXTERN PetscErrorCode SNESGetJacobian(SNES, Mat *, Mat *, SNESJacobianFn **, void **);
533: PETSC_EXTERN SNESFunctionFn SNESObjectiveComputeFunctionDefaultFD;
534: PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefault;
535: PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefaultColor;
536: PETSC_EXTERN PetscErrorCode SNESPruneJacobianColor(SNES, Mat, Mat);
537: PETSC_EXTERN PetscErrorCode SNESSetComputeInitialGuess(SNES, SNESInitialGuessFn *, void *);
538: PETSC_EXTERN PetscErrorCode SNESSetPicard(SNES, Vec, SNESFunctionFn *, Mat, Mat, SNESJacobianFn *, void *);
539: PETSC_EXTERN PetscErrorCode SNESGetPicard(SNES, Vec *, SNESFunctionFn **, Mat *, Mat *, SNESJacobianFn **, void **);
540: PETSC_EXTERN SNESFunctionFn SNESPicardComputeFunction;
541: PETSC_EXTERN SNESFunctionFn SNESPicardComputeMFFunction;
542: PETSC_EXTERN SNESJacobianFn SNESPicardComputeJacobian;

544: PETSC_EXTERN PetscErrorCode SNESSetObjective(SNES, SNESObjectiveFn *, void *);
545: PETSC_EXTERN PetscErrorCode SNESGetObjective(SNES, SNESObjectiveFn **, void **);
546: PETSC_EXTERN PetscErrorCode SNESComputeObjective(SNES, Vec, PetscReal *);

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

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

553:    Values:
554: +   `SNES_NORM_DEFAULT`            - use the default behavior for the current `SNESType`
555: .   `SNES_NORM_NONE`               - avoid all norm computations
556: .   `SNES_NORM_ALWAYS`             - compute the norms whenever possible
557: .   `SNES_NORM_INITIAL_ONLY`       - compute the norm only when the algorithm starts
558: .   `SNES_NORM_FINAL_ONLY`         - compute the norm only when the algorithm finishes
559: -   `SNES_NORM_INITIAL_FINAL_ONLY` - compute the norm at the start and end of the algorithm

561:    Level: advanced

563:    Notes:
564:    Support for these is highly dependent on the solver.

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

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

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

573: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
574:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
575: E*/
576: typedef enum {
577:   SNES_NORM_DEFAULT            = -1,
578:   SNES_NORM_NONE               = 0,
579:   SNES_NORM_ALWAYS             = 1,
580:   SNES_NORM_INITIAL_ONLY       = 2,
581:   SNES_NORM_FINAL_ONLY         = 3,
582:   SNES_NORM_INITIAL_FINAL_ONLY = 4
583: } SNESNormSchedule;
584: PETSC_EXTERN const char *const *const SNESNormSchedules;

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

589:    Level: advanced

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

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

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

600:    Level: advanced

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

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

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

611:    Level: advanced

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

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

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

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

626:    Level: advanced

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

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

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

640:    Level: advanced

642:    Note:
643:    This method combines the benefits of `SNES_NORM_INITIAL_ONLY` and `SNES_NORM_FINAL_ONLY`.

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

648: PETSC_EXTERN PetscErrorCode SNESSetNormSchedule(SNES, SNESNormSchedule);
649: PETSC_EXTERN PetscErrorCode SNESGetNormSchedule(SNES, SNESNormSchedule *);
650: PETSC_EXTERN PetscErrorCode SNESSetFunctionNorm(SNES, PetscReal);
651: PETSC_EXTERN PetscErrorCode SNESGetFunctionNorm(SNES, PetscReal *);
652: PETSC_EXTERN PetscErrorCode SNESGetUpdateNorm(SNES, PetscReal *);
653: PETSC_EXTERN PetscErrorCode SNESGetSolutionNorm(SNES, PetscReal *);

655: /*E
656:    SNESFunctionType - Type of function computed

658:    Values:
659: +  `SNES_FUNCTION_DEFAULT`          - the default behavior for the current `SNESType`
660: .  `SNES_FUNCTION_UNPRECONDITIONED` - the original function provided
661: -  `SNES_FUNCTION_PRECONDITIONED`   - the modification of the function by the preconditioner

663:    Level: advanced

665:    Note:
666:    Support for these is dependent on the solver.

668: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
669:           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
670: E*/
671: typedef enum {
672:   SNES_FUNCTION_DEFAULT          = -1,
673:   SNES_FUNCTION_UNPRECONDITIONED = 0,
674:   SNES_FUNCTION_PRECONDITIONED   = 1
675: } SNESFunctionType;
676: PETSC_EXTERN const char *const *const SNESFunctionTypes;

678: PETSC_EXTERN PetscErrorCode SNESSetFunctionType(SNES, SNESFunctionType);
679: PETSC_EXTERN PetscErrorCode SNESGetFunctionType(SNES, SNESFunctionType *);

681: PETSC_EXTERN PetscErrorCode SNESSetNGS(SNES, SNESNGSFn *, void *);
682: PETSC_EXTERN PetscErrorCode SNESGetNGS(SNES, SNESNGSFn **, void **);
683: PETSC_EXTERN PetscErrorCode SNESComputeNGS(SNES, Vec, Vec);

685: PETSC_EXTERN PetscErrorCode SNESNGSSetSweeps(SNES, PetscInt);
686: PETSC_EXTERN PetscErrorCode SNESNGSGetSweeps(SNES, PetscInt *);
687: PETSC_EXTERN PetscErrorCode SNESNGSSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt);
688: PETSC_EXTERN PetscErrorCode SNESNGSGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *);

690: PETSC_EXTERN PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES, PetscBool);
691: PETSC_EXTERN PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES, PetscBool *);

693: PETSC_EXTERN PetscErrorCode SNESShellGetContext(SNES, void *);
694: PETSC_EXTERN PetscErrorCode SNESShellSetContext(SNES, void *);
695: PETSC_EXTERN PetscErrorCode SNESShellSetSolve(SNES, PetscErrorCode (*)(SNES, Vec));

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

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

702:    Level: beginner

704: .seealso: [](ch_snes), `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNES`
705: S*/
706: typedef struct _p_LineSearch *SNESLineSearch;

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

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

721:    Level: beginner

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

727: .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetType()`, `SNES`
728: J*/
729: typedef const char *SNESLineSearchType;
730: #define SNESLINESEARCHBT        "bt"
731: #define SNESLINESEARCHNLEQERR   "nleqerr"
732: #define SNESLINESEARCHBASIC     "basic"
733: #define SNESLINESEARCHNONE      "none"
734: #define SNESLINESEARCHSECANT    "secant"
735: #define SNESLINESEARCHL2        PETSC_DEPRECATED_MACRO(3, 24, 0, "SNESLINESEARCHSECANT", ) "secant"
736: #define SNESLINESEARCHCP        "cp"
737: #define SNESLINESEARCHSHELL     "shell"
738: #define SNESLINESEARCHNCGLINEAR "ncglinear"
739: #define SNESLINESEARCHBISECTION "bisection"

741: PETSC_EXTERN PetscFunctionList SNESList;
742: PETSC_EXTERN PetscClassId      SNESLINESEARCH_CLASSID;
743: PETSC_EXTERN PetscFunctionList SNESLineSearchList;

745: #define SNES_LINESEARCH_ORDER_LINEAR    1
746: #define SNES_LINESEARCH_ORDER_QUADRATIC 2
747: #define SNES_LINESEARCH_ORDER_CUBIC     3

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

752:   Calling Sequence:
753: + snes  - `SNES` context
754: - u     - the vector to project to the bounds

756:   Level: advanced

758:   Note:
759:   The deprecated `SNESLineSearchVIProjectFunc` still works as a replacement for `SNESLineSearchVIProjectFn` *.

761: .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
762: S*/
763: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode             SNESLineSearchVIProjectFn(SNES snes, Vec u);
764: PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVIProjectFn *SNESLineSearchVIProjectFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVIProjectFn*", );

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

770:   Calling Sequence:
771: + snes  - `SNES` context
772: . f     - the vector to compute the norm of
773: . u     - the current solution, entries that are on the VI bounds are ignored
774: - fnorm - the resulting norm

776:   Level: advanced

778:   Note:
779:   The deprecated `SNESLineSearchVINormFunc` still works as a replacement for `SNESLineSearchVINormFn` *.

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

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

789:   Calling Sequence:
790: + snes  - `SNES` context
791: . f     - the function vector to compute the directional derivative with
792: . u     - the current solution, entries that are on the VI bounds are ignored
793: . y     - the direction to compute the directional derivative
794: - fty   - the resulting directional derivative

796:   Level: advanced

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

802: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchApplyFn(SNESLineSearch);
803: PETSC_EXTERN_TYPEDEF typedef SNESLineSearchApplyFn      *SNESLineSearchApplyFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
804: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchShellApplyFn(SNESLineSearch, void *);
805: PETSC_EXTERN_TYPEDEF typedef SNESLineSearchShellApplyFn *SNESLineSearchUserFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );

807: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate(MPI_Comm, SNESLineSearch *);
808: PETSC_EXTERN PetscErrorCode SNESLineSearchReset(SNESLineSearch);
809: PETSC_EXTERN PetscErrorCode SNESLineSearchView(SNESLineSearch, PetscViewer);
810: PETSC_EXTERN PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *);
811: PETSC_EXTERN PetscErrorCode SNESLineSearchGetType(SNESLineSearch, SNESLineSearchType *);
812: PETSC_EXTERN PetscErrorCode SNESLineSearchSetType(SNESLineSearch, SNESLineSearchType);
813: PETSC_EXTERN PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch);
814: PETSC_EXTERN PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch, PetscErrorCode (*)(SNES, Vec, Vec));
815: PETSC_EXTERN PetscErrorCode SNESLineSearchSetUp(SNESLineSearch);
816: PETSC_EXTERN PetscErrorCode SNESLineSearchApply(SNESLineSearch, Vec, Vec, PetscReal *, Vec);
817: PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch, Vec, Vec, PetscBool *);
818: PETSC_EXTERN PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *);
819: PETSC_EXTERN PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch, PetscInt);

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

823: PETSC_EXTERN PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx);
824: PETSC_EXTERN PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);

826: PETSC_EXTERN PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx);
827: PETSC_EXTERN PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);

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

831: PETSC_EXTERN PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn *, SNESLineSearchVINormFn *, SNESLineSearchVIDirDerivFn *);
832: PETSC_EXTERN PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn **, SNESLineSearchVINormFn **, SNESLineSearchVIDirDerivFn **);

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

838: /* set and get the parameters and vectors */
839: PETSC_EXTERN PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
840: PETSC_EXTERN PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscInt);

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

844: PETSC_EXTERN PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch, PetscReal *);
845: PETSC_EXTERN PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch, PetscReal);

847: PETSC_EXTERN PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch, PetscReal *);
848: PETSC_EXTERN PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch, PetscReal);

850: PETSC_EXTERN PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch, PetscInt *);
851: PETSC_EXTERN PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch, PetscInt);

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

856:   Values:
857: +  `SNES_LINESEARCH_SUCCEEDED`       - the line search succeeded
858: .  `SNES_LINESEARCH_FAILED_NANORINF` - a not a number of infinity appeared in the computions
859: .  `SNES_LINESEARCH_FAILED_DOMAIN`   - the function was evaluated outside of its domain, see `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
860: .  `SNES_LINESEARCH_FAILED_REDUCT`   - the linear search failed to get the requested decrease in its norm or objective
861: .  `SNES_LINESEARCH_FAILED_USER`     - used by `SNESLINESEARCHNLEQERR` to indicate the user changed the search direction inappropriately
862: -  `SNES_LINESEARCH_FAILED_FUNCTION` - indicates the maximum number of function evaluations allowed has been surpassed, `SNESConvergedReason` is also
863:                                        set to `SNES_DIVERGED_FUNCTION_COUNT`

865:    Level: intermediate

867:    Developer Note:
868:    Some of these reasons overlap with values of `SNESConvergedReason`

870: .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`,
871:           `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
872: E*/
873: typedef enum {
874:   SNES_LINESEARCH_SUCCEEDED,
875:   SNES_LINESEARCH_FAILED_NANORINF,
876:   SNES_LINESEARCH_FAILED_DOMAIN,
877:   SNES_LINESEARCH_FAILED_REDUCT, /* INSUFFICIENT REDUCTION */
878:   SNES_LINESEARCH_FAILED_USER,
879:   SNES_LINESEARCH_FAILED_FUNCTION
880: } SNESLineSearchReason;

882: PETSC_EXTERN PetscErrorCode SNESLineSearchGetReason(SNESLineSearch, SNESLineSearchReason *);
883: PETSC_EXTERN PetscErrorCode SNESLineSearchSetReason(SNESLineSearch, SNESLineSearchReason);

885: PETSC_EXTERN PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch, Vec *, Vec *, Vec *, Vec *, Vec *);
886: PETSC_EXTERN PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch, Vec, Vec, Vec, Vec, Vec);

888: PETSC_EXTERN PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *);
889: PETSC_EXTERN PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch, PetscReal, PetscReal, PetscReal);
890: PETSC_EXTERN PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch);
891: PETSC_EXTERN PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch, PetscBool);

893: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitor(SNESLineSearch);
894: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, void *), void *, PetscCtxDestroyFn *);
895: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch, const char[], const char[], const char[], PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *));
896: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch);
897: PETSC_EXTERN PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch, PetscViewer);
898: PETSC_EXTERN PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch, PetscViewer *);
899: PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch, PetscViewerAndFormat *);

901: PETSC_EXTERN PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch, const char[]);
902: PETSC_EXTERN PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch, const char *[]);

904: /* Shell interface functions */
905: PETSC_EXTERN PetscErrorCode SNESLineSearchShellSetApply(SNESLineSearch, SNESLineSearchShellApplyFn *, void *);
906: PETSC_EXTERN PetscErrorCode SNESLineSearchShellGetApply(SNESLineSearch, SNESLineSearchShellApplyFn **, void **);

908: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellSetApply()", ) static inline PetscErrorCode SNESLineSearchShellSetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn *f, void *ctx)
909: {
910:   return SNESLineSearchShellSetApply(ls, f, ctx);
911: }

913: PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellGetApply()", ) static inline PetscErrorCode SNESLineSearchShellGetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn **f, void **ctx)
914: {
915:   return SNESLineSearchShellGetApply(ls, f, ctx);
916: }

918: /* BT interface functions */
919: PETSC_EXTERN PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch, PetscReal);
920: PETSC_EXTERN PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch, PetscReal *);

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

925: /* Routines for VI solver */
926: PETSC_EXTERN PetscErrorCode SNESVISetVariableBounds(SNES, Vec, Vec);
927: PETSC_EXTERN PetscErrorCode SNESVIGetVariableBounds(SNES, Vec *, Vec *);
928: PETSC_EXTERN PetscErrorCode SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
929: PETSC_EXTERN PetscErrorCode SNESVIGetInactiveSet(SNES, IS *);
930: PETSC_EXTERN PetscErrorCode SNESVIGetActiveSetIS(SNES, Vec, Vec, IS *);
931: PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES, Vec, Vec, PetscReal *);
932: PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFtY(SNES, Vec, Vec, Vec, PetscScalar *);
933: PETSC_EXTERN PetscErrorCode SNESVISetRedundancyCheck(SNES, PetscErrorCode (*)(SNES, IS, IS *, void *), void *);
934: PETSC_EXTERN PetscErrorCode SNESVIComputeMeritFunction(Vec, PetscReal *, PetscReal *);
935: PETSC_EXTERN PetscErrorCode SNESVIComputeFunction(SNES, Vec, Vec, void *);
936: PETSC_EXTERN PetscErrorCode DMSetVI(DM, IS);
937: PETSC_EXTERN PetscErrorCode DMDestroyVI(DM);

939: PETSC_EXTERN PetscErrorCode SNESTestLocalMin(SNES);

941: /* Should this routine be private? */
942: PETSC_EXTERN PetscErrorCode SNESComputeJacobian(SNES, Vec, Mat, Mat);
943: PETSC_EXTERN PetscErrorCode SNESTestJacobian(SNES, PetscReal *, PetscReal *);
944: PETSC_EXTERN PetscErrorCode SNESTestFunction(SNES);

946: PETSC_EXTERN PetscErrorCode SNESSetDM(SNES, DM);
947: PETSC_EXTERN PetscErrorCode SNESGetDM(SNES, DM *);
948: PETSC_EXTERN PetscErrorCode SNESSetNPC(SNES, SNES);
949: PETSC_EXTERN PetscErrorCode SNESGetNPC(SNES, SNES *);
950: PETSC_EXTERN PetscErrorCode SNESHasNPC(SNES, PetscBool *);
951: PETSC_EXTERN PetscErrorCode SNESApplyNPC(SNES, Vec, Vec, Vec);
952: PETSC_EXTERN PetscErrorCode SNESGetNPCFunction(SNES, Vec, PetscReal *);
953: PETSC_EXTERN PetscErrorCode SNESComputeFunctionDefaultNPC(SNES, Vec, Vec);
954: PETSC_EXTERN PetscErrorCode SNESSetNPCSide(SNES, PCSide);
955: PETSC_EXTERN PetscErrorCode SNESGetNPCSide(SNES, PCSide *);
956: PETSC_EXTERN PetscErrorCode SNESSetLineSearch(SNES, SNESLineSearch);
957: PETSC_EXTERN PetscErrorCode SNESGetLineSearch(SNES, SNESLineSearch *);

959: PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESGetLineSearch()", ) static inline PetscErrorCode SNESGetSNESLineSearch(SNES snes, SNESLineSearch *ls)
960: {
961:   return SNESGetLineSearch(snes, ls);
962: }
963: PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESSetLineSearch()", ) static inline PetscErrorCode SNESSetSNESLineSearch(SNES snes, SNESLineSearch ls)
964: {
965:   return SNESSetLineSearch(snes, ls);
966: }

968: PETSC_EXTERN PetscErrorCode SNESSetUpMatrices(SNES);
969: PETSC_EXTERN PetscErrorCode DMSNESSetFunction(DM, SNESFunctionFn *, void *);
970: PETSC_EXTERN PetscErrorCode DMSNESGetFunction(DM, SNESFunctionFn **, void **);
971: PETSC_EXTERN PetscErrorCode DMSNESSetFunctionContextDestroy(DM, PetscCtxDestroyFn *);
972: PETSC_EXTERN PetscErrorCode DMSNESSetMFFunction(DM, SNESFunctionFn *, void *);
973: PETSC_EXTERN PetscErrorCode DMSNESSetNGS(DM, SNESNGSFn *, void *);
974: PETSC_EXTERN PetscErrorCode DMSNESGetNGS(DM, SNESNGSFn **, void **);
975: PETSC_EXTERN PetscErrorCode DMSNESSetJacobian(DM, SNESJacobianFn *, void *);
976: PETSC_EXTERN PetscErrorCode DMSNESGetJacobian(DM, SNESJacobianFn **, void **);
977: PETSC_EXTERN PetscErrorCode DMSNESSetJacobianContextDestroy(DM, PetscCtxDestroyFn *);
978: PETSC_EXTERN PetscErrorCode DMSNESSetPicard(DM, SNESFunctionFn *, SNESJacobianFn *, void *);
979: PETSC_EXTERN PetscErrorCode DMSNESGetPicard(DM, SNESFunctionFn **, SNESJacobianFn **, void **);
980: PETSC_EXTERN PetscErrorCode DMSNESSetObjective(DM, SNESObjectiveFn *, void *);
981: PETSC_EXTERN PetscErrorCode DMSNESGetObjective(DM, SNESObjectiveFn **, void **);
982: PETSC_EXTERN PetscErrorCode DMCopyDMSNES(DM, DM);

984: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionFn(DMDALocalInfo *, void *, void *, void *);
985: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianFn(DMDALocalInfo *, void *, Mat, Mat, void *);
986: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveFn(DMDALocalInfo *, void *, PetscReal *, void *);

988: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionVecFn(DMDALocalInfo *, Vec, Vec, void *);
989: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianVecFn(DMDALocalInfo *, Vec, Mat, Mat, void *);
990: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveVecFn(DMDALocalInfo *, Vec, PetscReal *, void *);

992: PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocal(DM, InsertMode, DMDASNESFunctionFn *, void *);
993: PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocal(DM, DMDASNESJacobianFn *, void *);
994: PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocal(DM, DMDASNESObjectiveFn *, void *);
995: PETSC_EXTERN PetscErrorCode DMDASNESSetPicardLocal(DM, InsertMode, DMDASNESFunctionFn *, DMDASNESJacobianFn, void *);

997: PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocalVec(DM, InsertMode, DMDASNESFunctionVecFn *, void *);
998: PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocalVec(DM, DMDASNESJacobianVecFn *, void *);
999: PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocalVec(DM, DMDASNESObjectiveVecFn *, void *);

1001: PETSC_EXTERN PetscErrorCode DMSNESSetBoundaryLocal(DM, PetscErrorCode (*)(DM, Vec, void *), void *);
1002: PETSC_EXTERN PetscErrorCode DMSNESSetObjectiveLocal(DM, PetscErrorCode (*)(DM, Vec, PetscReal *, void *), void *);
1003: PETSC_EXTERN PetscErrorCode DMSNESSetFunctionLocal(DM, PetscErrorCode (*)(DM, Vec, Vec, void *), void *);
1004: PETSC_EXTERN PetscErrorCode DMSNESSetJacobianLocal(DM, PetscErrorCode (*)(DM, Vec, Mat, Mat, void *), void *);
1005: PETSC_EXTERN PetscErrorCode DMSNESGetBoundaryLocal(DM, PetscErrorCode (**)(DM, Vec, void *), void **);
1006: PETSC_EXTERN PetscErrorCode DMSNESGetObjectiveLocal(DM, PetscErrorCode (**)(DM, Vec, PetscReal *, void *), void **);
1007: PETSC_EXTERN PetscErrorCode DMSNESGetFunctionLocal(DM, PetscErrorCode (**)(DM, Vec, Vec, void *), void **);
1008: PETSC_EXTERN PetscErrorCode DMSNESGetJacobianLocal(DM, PetscErrorCode (**)(DM, Vec, Mat, Mat, void *), void **);

1010: /* Routines for Multiblock solver */
1011: PETSC_EXTERN PetscErrorCode SNESMultiblockSetFields(SNES, const char[], PetscInt, const PetscInt *);
1012: PETSC_EXTERN PetscErrorCode SNESMultiblockSetIS(SNES, const char[], IS);
1013: PETSC_EXTERN PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
1014: PETSC_EXTERN PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);
1015: PETSC_EXTERN PetscErrorCode SNESMultiblockGetSubSNES(SNES, PetscInt *, SNES *[]);

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

1020:    Level: intermediate

1022: .seealso: [](ch_snes), `SNESMS`, `SNESMSGetType()`, `SNESMSSetType()`, `SNES`
1023: J*/
1024: typedef const char *SNESMSType;
1025: #define SNESMSM62       "m62"
1026: #define SNESMSEULER     "euler"
1027: #define SNESMSJAMESON83 "jameson83"
1028: #define SNESMSVLTP11    "vltp11"
1029: #define SNESMSVLTP21    "vltp21"
1030: #define SNESMSVLTP31    "vltp31"
1031: #define SNESMSVLTP41    "vltp41"
1032: #define SNESMSVLTP51    "vltp51"
1033: #define SNESMSVLTP61    "vltp61"

1035: PETSC_EXTERN PetscErrorCode SNESMSRegister(SNESMSType, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[]);
1036: PETSC_EXTERN PetscErrorCode SNESMSRegisterAll(void);
1037: PETSC_EXTERN PetscErrorCode SNESMSGetType(SNES, SNESMSType *);
1038: PETSC_EXTERN PetscErrorCode SNESMSSetType(SNES, SNESMSType);
1039: PETSC_EXTERN PetscErrorCode SNESMSGetDamping(SNES, PetscReal *);
1040: PETSC_EXTERN PetscErrorCode SNESMSSetDamping(SNES, PetscReal);
1041: PETSC_EXTERN PetscErrorCode SNESMSFinalizePackage(void);
1042: PETSC_EXTERN PetscErrorCode SNESMSInitializePackage(void);
1043: PETSC_EXTERN PetscErrorCode SNESMSRegisterDestroy(void);

1045: /*MC
1046:    SNESNGMRESRestartType - the restart approach used by `SNESNGMRES`

1048:   Values:
1049: +   `SNES_NGMRES_RESTART_NONE`       - never restart
1050: .   `SNES_NGMRES_RESTART_DIFFERENCE` - restart based upon difference criteria
1051: -   `SNES_NGMRES_RESTART_PERIODIC`   - restart after a fixed number of iterations

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

1057:    Level: intermediate

1059: .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1060:           `SNESNGMRESGetRestartType()`, `SNESNGMRESSelectType`
1061: M*/
1062: typedef enum {
1063:   SNES_NGMRES_RESTART_NONE       = 0,
1064:   SNES_NGMRES_RESTART_PERIODIC   = 1,
1065:   SNES_NGMRES_RESTART_DIFFERENCE = 2
1066: } SNESNGMRESRestartType;
1067: PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];

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

1073:    Values:
1074: +   `SNES_NGMRES_SELECT_NONE`       - choose the combined solution all the time
1075: .   `SNES_NGMRES_SELECT_DIFFERENCE` - choose based upon the selection criteria
1076: -   `SNES_NGMRES_SELECT_LINESEARCH` - choose based upon line search combination

1078:   Options Database Key:
1079: . -snes_ngmres_select_type<difference,none,linesearch> - select type

1081:    Level: intermediate

1083: .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1084:           `SNESNGMRESGetRestartType()`, `SNESNGMRESRestartType`
1085: M*/
1086: typedef enum {
1087:   SNES_NGMRES_SELECT_NONE       = 0,
1088:   SNES_NGMRES_SELECT_DIFFERENCE = 1,
1089:   SNES_NGMRES_SELECT_LINESEARCH = 2
1090: } SNESNGMRESSelectType;
1091: PETSC_EXTERN const char *const SNESNGMRESSelectTypes[];

1093: PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartType(SNES, SNESNGMRESRestartType);
1094: PETSC_EXTERN PetscErrorCode SNESNGMRESSetSelectType(SNES, SNESNGMRESSelectType);
1095: PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartFmRise(SNES, PetscBool);
1096: PETSC_EXTERN PetscErrorCode SNESNGMRESGetRestartFmRise(SNES, PetscBool *);

1098: /*MC
1099:    SNESNCGType - the conjugate update approach for `SNESNCG`

1101:    Values:
1102: +   `SNES_NCG_FR`  - Fletcher-Reeves update
1103: .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update, the default and the only one that tolerates generalized search directions
1104: .   `SNES_NCG_HS`  - Hestenes-Steifel update
1105: .   `SNES_NCG_DY`  - Dai-Yuan update
1106: -   `SNES_NCG_CD`  - Conjugate Descent update

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

1111:    Level: intermediate

1113: .seealso: `SNES, `SNESNCG`, `SNESNCGSetType()`
1114: M*/
1115: typedef enum {
1116:   SNES_NCG_FR  = 0,
1117:   SNES_NCG_PRP = 1,
1118:   SNES_NCG_HS  = 2,
1119:   SNES_NCG_DY  = 3,
1120:   SNES_NCG_CD  = 4
1121: } SNESNCGType;
1122: PETSC_EXTERN const char *const SNESNCGTypes[];

1124: PETSC_EXTERN PetscErrorCode SNESNCGSetType(SNES, SNESNCGType);

1126: /*MC
1127:    SNESQNScaleType - the scaling type used by `SNESQN`

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

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

1139:    Level: intermediate

1141: .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNRestartType`
1142: M*/
1143: typedef enum {
1144:   SNES_QN_SCALE_DEFAULT  = 0,
1145:   SNES_QN_SCALE_NONE     = 1,
1146:   SNES_QN_SCALE_SCALAR   = 2,
1147:   SNES_QN_SCALE_DIAGONAL = 3,
1148:   SNES_QN_SCALE_JACOBIAN = 4
1149: } SNESQNScaleType;
1150: PETSC_EXTERN const char *const SNESQNScaleTypes[];

1152: /*MC
1153:    SNESQNRestartType - the restart approached used by `SNESQN`

1155:    Values:
1156: +   `SNES_QN_RESTART_NONE`     - never restart
1157: .   `SNES_QN_RESTART_POWELL`   - restart based upon descent criteria
1158: -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations

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

1164:    Level: intermediate

1166: .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNScaleType`
1167: M*/
1168: typedef enum {
1169:   SNES_QN_RESTART_DEFAULT  = 0,
1170:   SNES_QN_RESTART_NONE     = 1,
1171:   SNES_QN_RESTART_POWELL   = 2,
1172:   SNES_QN_RESTART_PERIODIC = 3
1173: } SNESQNRestartType;
1174: PETSC_EXTERN const char *const SNESQNRestartTypes[];

1176: /*MC
1177:    SNESQNType - the type used by `SNESQN`

1179:   Values:
1180: +   `SNES_QN_LBFGS`      - LBFGS variant
1181: .   `SNES_QN_BROYDEN`    - Broyden variant
1182: -   `SNES_QN_BADBROYDEN` - Bad Broyden variant

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

1187:    Level: intermediate

1189: .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNSetType()`, `SNESQNScaleType`, `SNESQNRestartType`, `SNESQNSetRestartType()`
1190: M*/
1191: typedef enum {
1192:   SNES_QN_LBFGS      = 0,
1193:   SNES_QN_BROYDEN    = 1,
1194:   SNES_QN_BADBROYDEN = 2
1195: } SNESQNType;
1196: PETSC_EXTERN const char *const SNESQNTypes[];

1198: PETSC_EXTERN PetscErrorCode SNESQNSetType(SNES, SNESQNType);
1199: PETSC_EXTERN PetscErrorCode SNESQNSetScaleType(SNES, SNESQNScaleType);
1200: PETSC_EXTERN PetscErrorCode SNESQNSetRestartType(SNES, SNESQNRestartType);

1202: PETSC_EXTERN PetscErrorCode SNESNASMGetType(SNES, PCASMType *);
1203: PETSC_EXTERN PetscErrorCode SNESNASMSetType(SNES, PCASMType);
1204: PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomains(SNES, PetscInt *, SNES *[], VecScatter *[], VecScatter *[], VecScatter *[]);
1205: PETSC_EXTERN PetscErrorCode SNESNASMSetSubdomains(SNES, PetscInt, SNES[], VecScatter[], VecScatter[], VecScatter[]);
1206: PETSC_EXTERN PetscErrorCode SNESNASMSetDamping(SNES, PetscReal);
1207: PETSC_EXTERN PetscErrorCode SNESNASMGetDamping(SNES, PetscReal *);
1208: PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomainVecs(SNES, PetscInt *, Vec *[], Vec *[], Vec *[], Vec *[]);
1209: PETSC_EXTERN PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES, PetscBool);
1210: PETSC_EXTERN PetscErrorCode SNESNASMGetSNES(SNES, PetscInt, SNES *);
1211: PETSC_EXTERN PetscErrorCode SNESNASMGetNumber(SNES, PetscInt *);
1212: PETSC_EXTERN PetscErrorCode SNESNASMSetWeight(SNES, Vec);

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

1217:   Values:
1218: + `SNES_COMPOSITE_ADDITIVE`        - results from application of all preconditioners are added together
1219: . `SNES_COMPOSITE_MULTIPLICATIVE`  - preconditioners are applied sequentially to the residual freshly
1220:                                      computed after the previous preconditioner application
1221: - `SNES_COMPOSITE_ADDITIVEOPTIMAL` - uses a linear combination of the solutions obtained with each preconditioner that approximately minimize the function
1222:                                      value at the new iteration.

1224:    Level: beginner

1226: .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `PCCompositeType`
1227: E*/
1228: typedef enum {
1229:   SNES_COMPOSITE_ADDITIVE,
1230:   SNES_COMPOSITE_MULTIPLICATIVE,
1231:   SNES_COMPOSITE_ADDITIVEOPTIMAL
1232: } SNESCompositeType;
1233: PETSC_EXTERN const char *const SNESCompositeTypes[];

1235: PETSC_EXTERN PetscErrorCode SNESCompositeSetType(SNES, SNESCompositeType);
1236: PETSC_EXTERN PetscErrorCode SNESCompositeAddSNES(SNES, SNESType);
1237: PETSC_EXTERN PetscErrorCode SNESCompositeGetSNES(SNES, PetscInt, SNES *);
1238: PETSC_EXTERN PetscErrorCode SNESCompositeGetNumber(SNES, PetscInt *);
1239: PETSC_EXTERN PetscErrorCode SNESCompositeSetDamping(SNES, PetscInt, PetscReal);

1241: PETSC_EXTERN PetscErrorCode SNESPatchSetDiscretisationInfo(SNES, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
1242: PETSC_EXTERN PetscErrorCode SNESPatchSetComputeOperator(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
1243: PETSC_EXTERN PetscErrorCode SNESPatchSetComputeFunction(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
1244: PETSC_EXTERN PetscErrorCode SNESPatchSetConstructType(SNES, PCPatchConstructType, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *);
1245: PETSC_EXTERN PetscErrorCode SNESPatchSetCellNumbering(SNES, PetscSection);

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

1250:    Values:
1251: +  `SNES_FAS_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `SNESFASSetCycles()`
1252: .  `SNES_FAS_ADDITIVE`                 - additive FAS cycle
1253: .  `SNES_FAS_FULL`                     - full FAS cycle
1254: -  `SNES_FAS_KASKADE`                  - Kaskade FAS cycle

1256:    Level: beginner

1258: .seealso: [](ch_snes), `SNESFAS`, `PCMGSetType()`, `PCMGType`
1259: E*/
1260: typedef enum {
1261:   SNES_FAS_MULTIPLICATIVE,
1262:   SNES_FAS_ADDITIVE,
1263:   SNES_FAS_FULL,
1264:   SNES_FAS_KASKADE
1265: } SNESFASType;
1266: PETSC_EXTERN const char *const SNESFASTypes[];

1268: /* called on the finest level FAS instance*/
1269: PETSC_EXTERN PetscErrorCode SNESFASSetType(SNES, SNESFASType);
1270: PETSC_EXTERN PetscErrorCode SNESFASGetType(SNES, SNESFASType *);
1271: PETSC_EXTERN PetscErrorCode SNESFASSetLevels(SNES, PetscInt, MPI_Comm *);
1272: PETSC_EXTERN PetscErrorCode SNESFASGetLevels(SNES, PetscInt *);
1273: PETSC_EXTERN PetscErrorCode SNESFASGetCycleSNES(SNES, PetscInt, SNES *);
1274: PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothUp(SNES, PetscInt);
1275: PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothDown(SNES, PetscInt);
1276: PETSC_EXTERN PetscErrorCode SNESFASSetCycles(SNES, PetscInt);
1277: PETSC_EXTERN PetscErrorCode SNESFASSetMonitor(SNES, PetscViewerAndFormat *, PetscBool);
1278: PETSC_EXTERN PetscErrorCode SNESFASSetLog(SNES, PetscBool);

1280: PETSC_EXTERN PetscErrorCode SNESFASSetGalerkin(SNES, PetscBool);
1281: PETSC_EXTERN PetscErrorCode SNESFASGetGalerkin(SNES, PetscBool *);
1282: PETSC_EXTERN PetscErrorCode SNESFASGalerkinFunctionDefault(SNES, Vec, Vec, void *);

1284: /* called on any level -- "Cycle" FAS instance */
1285: PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmoother(SNES, SNES *);
1286: PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherUp(SNES, SNES *);
1287: PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherDown(SNES, SNES *);
1288: PETSC_EXTERN PetscErrorCode SNESFASCycleGetCorrection(SNES, SNES *);
1289: PETSC_EXTERN PetscErrorCode SNESFASCycleGetInterpolation(SNES, Mat *);
1290: PETSC_EXTERN PetscErrorCode SNESFASCycleGetRestriction(SNES, Mat *);
1291: PETSC_EXTERN PetscErrorCode SNESFASCycleGetInjection(SNES, Mat *);
1292: PETSC_EXTERN PetscErrorCode SNESFASCycleGetRScale(SNES, Vec *);
1293: PETSC_EXTERN PetscErrorCode SNESFASCycleSetCycles(SNES, PetscInt);
1294: PETSC_EXTERN PetscErrorCode SNESFASCycleIsFine(SNES, PetscBool *);

1296: /* called on the (outer) finest level FAS to set/get parameters on any level instance */
1297: PETSC_EXTERN PetscErrorCode SNESFASSetInterpolation(SNES, PetscInt, Mat);
1298: PETSC_EXTERN PetscErrorCode SNESFASGetInterpolation(SNES, PetscInt, Mat *);
1299: PETSC_EXTERN PetscErrorCode SNESFASSetRestriction(SNES, PetscInt, Mat);
1300: PETSC_EXTERN PetscErrorCode SNESFASGetRestriction(SNES, PetscInt, Mat *);
1301: PETSC_EXTERN PetscErrorCode SNESFASSetInjection(SNES, PetscInt, Mat);
1302: PETSC_EXTERN PetscErrorCode SNESFASGetInjection(SNES, PetscInt, Mat *);
1303: PETSC_EXTERN PetscErrorCode SNESFASSetRScale(SNES, PetscInt, Vec);
1304: PETSC_EXTERN PetscErrorCode SNESFASGetRScale(SNES, PetscInt, Vec *);
1305: PETSC_EXTERN PetscErrorCode SNESFASSetContinuation(SNES, PetscBool);

1307: PETSC_EXTERN PetscErrorCode SNESFASGetSmoother(SNES, PetscInt, SNES *);
1308: PETSC_EXTERN PetscErrorCode SNESFASGetSmootherUp(SNES, PetscInt, SNES *);
1309: PETSC_EXTERN PetscErrorCode SNESFASGetSmootherDown(SNES, PetscInt, SNES *);
1310: PETSC_EXTERN PetscErrorCode SNESFASGetCoarseSolve(SNES, SNES *);

1312: /* parameters for full FAS */
1313: PETSC_EXTERN PetscErrorCode SNESFASFullSetDownSweep(SNES, PetscBool);
1314: PETSC_EXTERN PetscErrorCode SNESFASCreateCoarseVec(SNES, Vec *);
1315: PETSC_EXTERN PetscErrorCode SNESFASRestrict(SNES, Vec, Vec);
1316: PETSC_EXTERN PetscErrorCode SNESFASFullSetTotal(SNES, PetscBool);
1317: PETSC_EXTERN PetscErrorCode SNESFASFullGetTotal(SNES, PetscBool *);

1319: PETSC_EXTERN PetscErrorCode DMPlexSetSNESVariableBounds(DM, SNES);
1320: PETSC_EXTERN PetscErrorCode DMSNESCheckDiscretization(SNES, DM, PetscReal, Vec, PetscReal, PetscReal[]);
1321: PETSC_EXTERN PetscErrorCode DMSNESCheckResidual(SNES, DM, Vec, PetscReal, PetscReal *);
1322: PETSC_EXTERN PetscErrorCode DMSNESCheckJacobian(SNES, DM, Vec, PetscReal, PetscBool *, PetscReal *);
1323: PETSC_EXTERN PetscErrorCode DMSNESCheckFromOptions(SNES, Vec);
1324: PETSC_EXTERN PetscErrorCode DMSNESComputeJacobianAction(DM, Vec, Vec, Vec, void *);
1325: PETSC_EXTERN PetscErrorCode DMSNESCreateJacobianMF(DM, Vec, void *, Mat *);

1327: PETSC_EXTERN PetscErrorCode SNESNewtonALSetFunction(SNES, SNESFunctionFn *, void *ctx);
1328: PETSC_EXTERN PetscErrorCode SNESNewtonALGetFunction(SNES, SNESFunctionFn **, void **ctx);
1329: PETSC_EXTERN PetscErrorCode SNESNewtonALComputeFunction(SNES, Vec, Vec);
1330: PETSC_EXTERN PetscErrorCode SNESNewtonALGetLoadParameter(SNES, PetscReal *);

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

1339:    Values:
1340: +   `SNES_NEWTONAL_CORRECTION_EXACT` - choose the correction which exactly satisfies the constraint
1341: -   `SNES_NEWTONAL_CORRECTION_NORMAL` - choose the correction in the updated normal hyper-surface to the constraint surface

1343:    Options Database Key:
1344: . -snes_newtonal_correction_type <exact> - select type from <exact,normal>

1346:    Level: intermediate

1348: .seealso: `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetCorrectionType()`
1349: M*/
1350: typedef enum {
1351:   SNES_NEWTONAL_CORRECTION_EXACT  = 0,
1352:   SNES_NEWTONAL_CORRECTION_NORMAL = 1,
1353: } SNESNewtonALCorrectionType;
1354: PETSC_EXTERN const char *const SNESNewtonALCorrectionTypes[];

1356: PETSC_EXTERN PetscErrorCode SNESNewtonALSetCorrectionType(SNES, SNESNewtonALCorrectionType);