xref: /petsc/include/petscsnes.h (revision 76c6338944e4871467ad7a763eee41b62845c2b0)
1f26ada1bSBarry Smith /*
2eef9c623SLois Curfman McInnes     User interface for the nonlinear solvers package.
3f26ada1bSBarry Smith */
4a4963045SJacob Faibussowitsch #pragma once
5ac09b921SBarry Smith 
62c8e378dSBarry Smith #include <petscksp.h>
71e25c274SJed Brown #include <petscdmtypes.h>
87ed8fce4SMatthew G. Knepley #include <petscfvtypes.h>
91e25c274SJed Brown #include <petscdmdatypes.h>
10e5148a0bSMatthew G. Knepley #include <petscsnestypes.h>
11b1f5cb9dSBarry Smith 
12ac09b921SBarry Smith /* SUBMANSEC = SNES */
13ac09b921SBarry Smith 
1476bdecfbSBarry Smith /*J
150b4b7b1cSBarry Smith    SNESType - String with the name of a PETSc `SNES` method. These are all the nonlinear solvers that PETSc provides.
16435da068SBarry Smith 
17435da068SBarry Smith    Level: beginner
18435da068SBarry Smith 
190b4b7b1cSBarry Smith    Note:
200b4b7b1cSBarry Smith    Use `SNESSetType()` or the options database key `-snes_type` to set the specific nonlinear solver algorithm to use with a given `SNES` object
210b4b7b1cSBarry Smith 
221cc06b55SBarry Smith .seealso: [](doc_nonlinsolve), [](ch_snes), `SNESSetType()`, `SNES`, `SNESCreate()`, `SNESDestroy()`, `SNESSetFromOptions()`
2376bdecfbSBarry Smith J*/
2419fd82e9SBarry Smith typedef const char *SNESType;
2504d7464bSBarry Smith #define SNESNEWTONLS         "newtonls"
2604d7464bSBarry Smith #define SNESNEWTONTR         "newtontr"
2741ba4c6cSHeeho Park #define SNESNEWTONTRDC       "newtontrdc"
281d6018f0SLisandro Dalcin #define SNESPYTHON           "python"
29d5c3842bSBarry Smith #define SNESNRICHARDSON      "nrichardson"
30b79b07cfSJed Brown #define SNESKSPONLY          "ksponly"
311ef27442SStefano Zampini #define SNESKSPTRANSPOSEONLY "ksptransposeonly"
32f450aa47SBarry Smith #define SNESVINEWTONRSLS     "vinewtonrsls"
33f450aa47SBarry Smith #define SNESVINEWTONSSLS     "vinewtonssls"
344a0c5b0cSMatthew G Knepley #define SNESNGMRES           "ngmres"
354b11644fSPeter Brune #define SNESQN               "qn"
36c5ae4b9aSBarry Smith #define SNESSHELL            "shell"
37be95d8f1SBarry Smith #define SNESNGS              "ngs"
38fef7b6d8SPeter Brune #define SNESNCG              "ncg"
39421d9b32SPeter Brune #define SNESFAS              "fas"
4037e1895aSJed Brown #define SNESMS               "ms"
41eaedb033SPeter Brune #define SNESNASM             "nasm"
42f31c9d25SPeter Brune #define SNESANDERSON         "anderson"
43d728fb7dSPeter Brune #define SNESASPIN            "aspin"
44eed5f15bSPeter Brune #define SNESCOMPOSITE        "composite"
45561742edSMatthew G. Knepley #define SNESPATCH            "patch"
4697276fddSZach Atkins #define SNESNEWTONAL         "newtonal"
47d5c3842bSBarry Smith 
48deeb6e72SMatthew Knepley /* Logging support */
49014dd563SJed Brown PETSC_EXTERN PetscClassId SNES_CLASSID;
5022c6f798SBarry Smith PETSC_EXTERN PetscClassId DMSNES_CLASSID;
51deeb6e72SMatthew Knepley 
52607a6623SBarry Smith PETSC_EXTERN PetscErrorCode SNESInitializePackage(void);
534bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESFinalizePackage(void);
54deeb6e72SMatthew Knepley 
55014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate(MPI_Comm, SNES *);
5677e5a1f9SBarry Smith PETSC_EXTERN PetscErrorCode SNESParametersInitialize(SNES);
57014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESReset(SNES);
58014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESDestroy(SNES *);
5919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetType(SNES, SNESType);
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMonitor(SNES, PetscInt, PetscReal);
6149abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorSet(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, void *), void *, PetscCtxDestroyFn *);
62d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorSetFromOptions(SNES, const char[], const char[], const char[], PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(SNES, PetscViewerAndFormat *));
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMonitorCancel(SNES);
648a70d858SHong Zhang PETSC_EXTERN PetscErrorCode SNESMonitorSAWs(SNES, PetscInt, PetscReal, void *);
658a70d858SHong Zhang PETSC_EXTERN PetscErrorCode SNESMonitorSAWsCreate(SNES, void **);
668a70d858SHong Zhang PETSC_EXTERN PetscErrorCode SNESMonitorSAWsDestroy(void **);
67014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetConvergenceHistory(SNES, PetscReal[], PetscInt[], PetscInt, PetscBool);
68014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetConvergenceHistory(SNES, PetscReal *[], PetscInt *[], PetscInt *);
69014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetUp(SNES);
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSolve(SNES, Vec, Vec);
71014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetErrorIfNotConverged(SNES, PetscBool);
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetErrorIfNotConverged(SNES, PetscBool *);
732d157150SStefano Zampini PETSC_EXTERN PetscErrorCode SNESConverged(SNES, PetscInt, PetscReal, PetscReal, PetscReal);
74e113a28aSBarry Smith 
75fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetWorkVecs(SNES, PetscInt);
7684cb2905SBarry Smith 
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*)(SNES));
784d8f6ca9SMatthew Knepley 
79bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode SNESRegister(const char[], PetscErrorCode (*)(SNES));
8030de9b25SBarry Smith 
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetKSP(SNES, KSP *);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetKSP(SNES, KSP);
833cd8a7caSMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetSolution(SNES, Vec);
84014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetSolution(SNES, Vec *);
85014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetSolutionUpdate(SNES, Vec *);
86014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetRhs(SNES, Vec *);
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESView(SNES, PetscViewer);
8855849f57SBarry Smith PETSC_EXTERN PetscErrorCode SNESLoad(SNES, PetscViewer);
8949abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode SNESConvergedReasonViewSet(SNES, PetscErrorCode (*)(SNES, void *), void *, PetscCtxDestroyFn *);
90fe2efc57SMark PETSC_EXTERN PetscErrorCode SNESViewFromOptions(SNES, PetscObject, const char[]);
9119a666eeSBarry Smith PETSC_EXTERN PetscErrorCode SNESConvergedReasonView(SNES, PetscViewer);
9219a666eeSBarry Smith PETSC_EXTERN PetscErrorCode SNESConvergedReasonViewFromOptions(SNES);
93c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode SNESConvergedReasonViewCancel(SNES);
9419a666eeSBarry Smith 
95edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SNESConvergedReasonView()", ) static inline PetscErrorCode SNESReasonView(SNES snes, PetscViewer v)
96d71ae5a4SJacob Faibussowitsch {
979371c9d4SSatish Balay   return SNESConvergedReasonView(snes, v);
989371c9d4SSatish Balay }
99edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "SNESConvergedReasonViewFromOptions()", ) static inline PetscErrorCode SNESReasonViewFromOptions(SNES snes)
100d71ae5a4SJacob Faibussowitsch {
1019371c9d4SSatish Balay   return SNESConvergedReasonViewFromOptions(snes);
1029371c9d4SSatish Balay }
10355849f57SBarry Smith 
10455849f57SBarry Smith #define SNES_FILE_CLASSID 1211224
1057bc3d0afSSatish Balay 
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetOptionsPrefix(SNES, const char[]);
107014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESAppendOptionsPrefix(SNES, const char[]);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetOptionsPrefix(SNES, const char *[]);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetFromOptions(SNES);
11087d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESResetFromOptions(SNES);
11140191667SLois Curfman McInnes 
1123565c898SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetUseMatrixFree(SNES, PetscBool, PetscBool);
1133565c898SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetUseMatrixFree(SNES, PetscBool *, PetscBool *);
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSNESMF(SNES, Mat *);
115bc13fc8dSBarry Smith PETSC_EXTERN PetscErrorCode MatSNESMFGetSNES(Mat, SNES *);
116208be567SBarry Smith PETSC_EXTERN PetscErrorCode MatSNESMFSetReuseBase(Mat, PetscBool);
117208be567SBarry Smith PETSC_EXTERN PetscErrorCode MatSNESMFGetReuseBase(Mat, PetscBool *);
118d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDComputeJacobian(SNES, Vec, Mat, Mat, void *);
119f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode MatCreateSNESMFMore(SNES, Vec, Mat *);
120f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode MatSNESMFMoreSetParameters(Mat, PetscReal, PetscReal, PetscReal);
1218f6e3e37SBarry Smith 
12219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetType(SNES, SNESType *);
123fbcc4530SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESMonitorDefaultSetUp(SNES, PetscViewerAndFormat *);
124d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorDefault(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
1251f60017eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorScaling(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
126d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorRange(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
127d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorRatio(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
128d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorRatioSetUp(SNES, PetscViewerAndFormat *);
129d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorSolution(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
130d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorResidual(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
131d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorSolutionUpdate(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
132d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorDefaultShort(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
133d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorDefaultField(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
134d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorJacUpdateSpectrum(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
135d43b4f6eSBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorFields(SNES, PetscInt, PetscReal, PetscViewerAndFormat *);
136798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSNESResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
137798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSNESResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *);
138798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSNESResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **);
139e5f7ee39SBarry Smith 
140014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt, PetscInt);
141e4d06f11SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESSetDivergenceTolerance(SNES, PetscReal);
142014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *, PetscInt *);
143e4d06f11SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESGetDivergenceTolerance(SNES, PetscReal *);
14485216dc7SFande Kong PETSC_EXTERN PetscErrorCode SNESGetForceIteration(SNES, PetscBool *);
145be5caee7SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetForceIteration(SNES, PetscBool);
146014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetIterationNumber(SNES, PetscInt *);
147014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetIterationNumber(SNES, PetscInt);
1483d4c4710SBarry Smith 
1494a221d59SStefano Zampini /*E
1504a221d59SStefano Zampini    SNESNewtonTRFallbackType - type of fallback in case the solution of the trust-region subproblem is outside of the radius
1514a221d59SStefano Zampini 
1524a221d59SStefano Zampini    Values:
1534a221d59SStefano Zampini +  `SNES_TR_FALLBACK_NEWTON` - use scaled Newton step
1544a221d59SStefano Zampini .  `SNES_TR_FALLBACK_CAUCHY` - use Cauchy direction
1554a221d59SStefano Zampini -  `SNES_TR_FALLBACK_DOGLEG` - use dogleg method
15616a05f60SBarry Smith 
15716a05f60SBarry Smith    Level: intermediate
15816a05f60SBarry Smith 
1591cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNEWTONTRDC`
1604a221d59SStefano Zampini E*/
1614a221d59SStefano Zampini typedef enum {
1624a221d59SStefano Zampini   SNES_TR_FALLBACK_NEWTON,
1634a221d59SStefano Zampini   SNES_TR_FALLBACK_CAUCHY,
1644a221d59SStefano Zampini   SNES_TR_FALLBACK_DOGLEG,
1654a221d59SStefano Zampini } SNESNewtonTRFallbackType;
1664a221d59SStefano Zampini 
1674a221d59SStefano Zampini PETSC_EXTERN const char *const SNESNewtonTRFallbackTypes[];
1684a221d59SStefano Zampini 
169c9368356SGlenn Hammond PETSC_EXTERN PetscErrorCode SNESNewtonTRSetPreCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, PetscBool *, void *), void *ctx);
170c9368356SGlenn Hammond PETSC_EXTERN PetscErrorCode SNESNewtonTRGetPreCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, PetscBool *, void *), void **ctx);
171c9368356SGlenn Hammond PETSC_EXTERN PetscErrorCode SNESNewtonTRSetPostCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);
172c9368356SGlenn Hammond PETSC_EXTERN PetscErrorCode SNESNewtonTRGetPostCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);
1734a221d59SStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRSetFallbackType(SNES, SNESNewtonTRFallbackType);
1744a221d59SStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRPreCheck(SNES, Vec, Vec, PetscBool *);
1754a221d59SStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRPostCheck(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *);
17624fb275aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRSetNormType(SNES, NormType);
17724fb275aSStefano Zampini 
17824fb275aSStefano Zampini /*E
17924fb275aSStefano Zampini     SNESNewtonTRQNType - type of quasi-Newton model to use
18024fb275aSStefano Zampini 
18124fb275aSStefano Zampini    Values:
18224fb275aSStefano Zampini +  `SNES_TR_QN_NONE` - do not use a quasi-Newton model
18324fb275aSStefano Zampini .  `SNES_TR_QN_SAME` - use the same quasi-Newton model for matrix and preconditioner
18424fb275aSStefano Zampini -  `SNES_TR_QN_DIFFERENT` - use different quasi-Newton models for matrix and preconditioner
18524fb275aSStefano Zampini 
18624fb275aSStefano Zampini    Level: intermediate
18724fb275aSStefano Zampini 
18824fb275aSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`
18924fb275aSStefano Zampini E*/
19024fb275aSStefano Zampini typedef enum {
19124fb275aSStefano Zampini   SNES_TR_QN_NONE,
19224fb275aSStefano Zampini   SNES_TR_QN_SAME,
19324fb275aSStefano Zampini   SNES_TR_QN_DIFFERENT,
19424fb275aSStefano Zampini } SNESNewtonTRQNType;
19524fb275aSStefano Zampini 
19624fb275aSStefano Zampini PETSC_EXTERN const char *const SNESNewtonTRQNTypes[];
19724fb275aSStefano Zampini 
19824fb275aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRSetQNType(SNES, SNESNewtonTRQNType);
1997cb011f5SBarry Smith 
2003201ab8dSStefano Zampini PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 22, 0, "SNESNewtonTRSetTolerances()", ) PetscErrorCode SNESSetTrustRegionTolerance(SNES, PetscReal);
2013201ab8dSStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRSetTolerances(SNES, PetscReal, PetscReal, PetscReal);
2023201ab8dSStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *);
2033201ab8dSStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRSetUpdateParameters(SNES, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
2043201ab8dSStefano Zampini PETSC_EXTERN PetscErrorCode SNESNewtonTRGetUpdateParameters(SNES, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
2053201ab8dSStefano Zampini 
20641ba4c6cSHeeho Park PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES, PetscBool *);
20741ba4c6cSHeeho Park PETSC_EXTERN PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, PetscBool *, void *), void *ctx);
20841ba4c6cSHeeho Park PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, PetscBool *, void *), void **ctx);
20941ba4c6cSHeeho Park PETSC_EXTERN PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);
21041ba4c6cSHeeho Park PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);
21141ba4c6cSHeeho Park 
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetNonlinearStepFailures(SNES, PetscInt *);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES, PetscInt);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES, PetscInt *);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetNumberFunctionEvals(SNES, PetscInt *);
216b850b91aSLisandro Dalcin 
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLagPreconditioner(SNES, PetscInt);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLagPreconditioner(SNES, PetscInt *);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLagJacobian(SNES, PetscInt);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLagJacobian(SNES, PetscInt *);
22137ec4e1aSPeter Brune PETSC_EXTERN PetscErrorCode SNESSetLagPreconditionerPersists(SNES, PetscBool);
22237ec4e1aSPeter Brune PETSC_EXTERN PetscErrorCode SNESSetLagJacobianPersists(SNES, PetscBool);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetGridSequence(SNES, PetscInt);
224fa19ca70SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetGridSequence(SNES, PetscInt *);
225a8054027SBarry Smith 
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLinearSolveIterations(SNES, PetscInt *);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLinearSolveFailures(SNES, PetscInt *);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetMaxLinearSolveFailures(SNES, PetscInt);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetMaxLinearSolveFailures(SNES, PetscInt *);
230971e163fSPeter Brune PETSC_EXTERN PetscErrorCode SNESSetCountersReset(SNES, PetscBool);
23112b1dd1aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESResetCounters(SNES);
2323d4c4710SBarry Smith 
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPSetUseEW(SNES, PetscBool);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPGetUseEW(SNES, PetscBool *);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPSetParametersEW(SNES, PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPGetParametersEW(SNES, PetscInt *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
237eafb4bcbSBarry Smith 
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMonitorLGRange(SNES, PetscInt, PetscReal, void *);
239eafb4bcbSBarry Smith 
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetApplicationContext(SNES, void *);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetApplicationContext(SNES, void *);
24249abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode SNESSetComputeApplicationContext(SNES, PetscErrorCode (*)(SNES, void **), PetscCtxDestroyFn *);
243184914b5SBarry Smith 
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESPythonSetType(SNES, const char[]);
245ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode SNESPythonGetType(SNES, const char *[]);
2461d6018f0SLisandro Dalcin 
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetFunctionDomainError(SNES);
248*76c63389SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetObjectiveDomainError(SNES);
249cc6b0f04SFande Kong PETSC_EXTERN PetscErrorCode SNESSetJacobianDomainError(SNES);
250b351a90bSFande Kong PETSC_EXTERN PetscErrorCode SNESSetCheckJacobianDomainError(SNES, PetscBool);
2518383d7d7SFande Kong PETSC_EXTERN PetscErrorCode SNESGetCheckJacobianDomainError(SNES, PetscBool *);
2526a388c36SPeter Brune 
253edd03b47SJacob Faibussowitsch #define SNES_CONVERGED_TR_DELTA_DEPRECATED SNES_CONVERGED_TR_DELTA PETSC_DEPRECATED_ENUM(3, 12, 0, "SNES_DIVERGED_TR_DELTA", )
254435da068SBarry Smith /*E
255dc4c0fb0SBarry Smith     SNESConvergedReason - reason a `SNESSolve()` was determined to have converged or diverged
256435da068SBarry Smith 
257dc4c0fb0SBarry Smith    Values:
258*76c63389SBarry Smith +  `SNES_CONVERGED_FNORM_ABS`         - $ ||F|| \le abstol $
259*76c63389SBarry Smith .  `SNES_CONVERGED_FNORM_RELATIVE`    - $ ||F|| <= rtol*||F(x_0))|| $ where $x_0 $ is the initial guess
260*76c63389SBarry Smith .  `SNES_CONVERGED_SNORM_RELATIVE`    - The 2-norm of the last step $ \le stol * ||x|| $ where $ x $ is the current solution
261379ef8d2SAlexander .  `SNES_CONVERGED_USER`              - The user has indicated convergence for an arbitrary reason
262dc4c0fb0SBarry Smith .  `SNES_DIVERGED_FUNCTION_COUNT`     - The user provided function has been called more times than the maximum set in `SNESSetTolerances()`
263dc4c0fb0SBarry Smith .  `SNES_DIVERGED_DTOL`               - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
264*76c63389SBarry Smith .  `SNES_DIVERGED_FUNCTION_NANORINF`  - the 2-norm of the current function evaluation is not-a-number (NaN) or infinity, (this
265*76c63389SBarry Smith                                         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*76c63389SBarry Smith .  `SNES_DIVERGED_OBJECTIVE_NANORINF` - the object function evaluation is not-a-number (NaN) or infinity, (this
267*76c63389SBarry Smith                                         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*76c63389SBarry Smith .  `SNES_DIVERGED_FUNCTION_DOMAIN`    - the function evaluation occurred outside the function's domain (function callback provided by
269*76c63389SBarry Smith                                         `SNESSetFunction()` called `SNESSetObjectiveDomainError()`) and the solver could not recover from this (by, for example, cutting the step size)
270*76c63389SBarry Smith .  `SNES_DIVERGED_OBJECTIVE_DOMAIN`   - the object function evaluation occurred outside the function's domain (function callback provided by
271*76c63389SBarry Smith                                         `SNESSetObjective()` called `SNESSetObjectiveDomainError()`) and the solver could not recover from this (by, for example, cutting the step size)
272*76c63389SBarry Smith .  `SNES_DIVERGED_JACOBIAN_DOMAIN`    - the Jacobian evaluation occurred outside the function's domain (function callback provided by
273*76c63389SBarry Smith                                         `SNESSetJacobian()` called `SNESSetJacobianDomainError()`)
274dc4c0fb0SBarry Smith .  `SNES_DIVERGED_MAX_IT`             - `SNESSolve()` has reached the maximum number of iterations requested
275dc4c0fb0SBarry Smith .  `SNES_DIVERGED_LINE_SEARCH`        - The line search has failed. This only occurs for `SNES` solvers that use a line search
276dc4c0fb0SBarry Smith .  `SNES_DIVERGED_LOCAL_MIN`          - the algorithm seems to have stagnated at a local minimum that is not zero.
277379ef8d2SAlexander -  `SNES_CONVERGED_ITERATING          - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
278f203c74bSBarry Smith 
27916a05f60SBarry Smith    Level: beginner
28016a05f60SBarry Smith 
281dc4c0fb0SBarry Smith     Notes:
2820b4b7b1cSBarry Smith    The two most common reasons for divergence are an incorrectly coded or computed Jacobian or failure or lack of convergence in the linear system
2830b4b7b1cSBarry Smith    (in this case we recommend
284dc4c0fb0SBarry Smith    testing with `-pc_type lu` to eliminate the linear solver as the cause of the problem).
285dc4c0fb0SBarry Smith 
2860b4b7b1cSBarry Smith    `SNES_DIVERGED_LOCAL_MIN` can only occur when using a `SNES` solver that uses a line search (`SNESLineSearch`).
287*76c63389SBarry Smith    The line search wants to $ \min Q(\alpha) = 1/2 || F(x + \alpha s) ||^2_2 $  this occurs
288*76c63389SBarry Smith    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*76c63389SBarry Smith    $ Q'(\alpha) = -F(x)^T F'(x)^(-1)^T F'(x+\alpha s)F(x+\alpha s)$; when $\alpha = 0$
290*76c63389SBarry Smith    $Q'(0) = - ||F(x)||^2_2 $ which is always NEGATIVE if $F'(x)$ is invertible. This means the Newton
291*76c63389SBarry Smith    direction is a descent direction and the line search should succeed if $\alpha $ is small enough.
292f203c74bSBarry Smith 
293*76c63389SBarry Smith    If $F'(x)$ is NOT invertible AND $F'(x)^T F(x) = 0 $ then $Q'(0) = 0 $ and the Newton direction
294f203c74bSBarry Smith    is NOT a descent direction so the line search will fail. All one can do at this point
295f203c74bSBarry Smith    is change the initial guess and try again.
296f203c74bSBarry Smith 
297f203c74bSBarry Smith    An alternative explanation: Newton's method can be regarded as replacing the function with
298*76c63389SBarry Smith    its linear approximation and minimizing the 2-norm of that. That is $F(x+s) \approx F(x) + F'(x)s$
299*76c63389SBarry Smith    so we minimize $ || F(x) + F'(x) s ||^2_2$ using Least Squares. If $F'(x)$ is invertible then
300*76c63389SBarry Smith    $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*76c63389SBarry Smith    exists a nontrivial (that is $F'(x)s \ne 0$) solution to the equation and this direction is
302*76c63389SBarry Smith    $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*76c63389SBarry Smith    = - (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*76c63389SBarry Smith    and $F'(x)^T F'(x)$ has no negative eigenvalues $Q'(0) < 0$ so $s$ is a descent direction and the line
305*76c63389SBarry Smith    search should succeed for small enough $\alpha$.
306f203c74bSBarry Smith 
307f203c74bSBarry Smith    Note that this RARELY happens in practice. Far more likely the linear system is not being solved
308f203c74bSBarry Smith    (well enough?) or the Jacobian is wrong.
309f203c74bSBarry Smith 
31087497f52SBarry Smith    `SNES_DIVERGED_MAX_IT` means that the solver reached the maximum number of iterations without satisfying any
31187497f52SBarry Smith    convergence criteria. `SNES_CONVERGED_ITS` means that `SNESConvergedSkip()` was chosen as the convergence test;
3124d0a8057SBarry Smith    thus the usual convergence criteria have not been checked and may or may not be satisfied.
313f203c74bSBarry Smith 
3141cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`, `SNESSetTolerances()`
315435da068SBarry Smith E*/
316184914b5SBarry Smith typedef enum {                       /* converged */
31701b82886SBarry Smith   SNES_CONVERGED_FNORM_ABS      = 2, /* ||F|| < atol */
31801b82886SBarry Smith   SNES_CONVERGED_FNORM_RELATIVE = 3, /* ||F|| < rtol*||F_initial|| */
3195358d0d4SBarry Smith   SNES_CONVERGED_SNORM_RELATIVE = 4, /* Newton computed step size small; || delta x || < stol || x || */
3203f149594SLisandro Dalcin   SNES_CONVERGED_ITS            = 5, /* maximum iterations reached */
321ce78bad3SBarry Smith   SNES_BREAKOUT_INNER_ITER      = 6, /* Flag to break out of inner loop after checking custom convergence, used in multi-phase flow when state changes */
322379ef8d2SAlexander   SNES_CONVERGED_USER           = 7, /* The user has indicated convergence for an arbitrary reason */
323184914b5SBarry Smith   /* diverged */
32446a9e3ceSBarry Smith   SNES_DIVERGED_FUNCTION_DOMAIN      = -1, /* the new x location passed the function is not in the domain of F */
325184914b5SBarry Smith   SNES_DIVERGED_FUNCTION_COUNT       = -2,
32646a9e3ceSBarry Smith   SNES_DIVERGED_LINEAR_SOLVE         = -3, /* the linear solve failed */
327*76c63389SBarry Smith   SNES_DIVERGED_FUNCTION_NANORINF    = -4,
328184914b5SBarry Smith   SNES_DIVERGED_MAX_IT               = -5,
329647a2e1fSBarry Smith   SNES_DIVERGED_LINE_SEARCH          = -6,  /* the line search failed */
3301e633543SBarry Smith   SNES_DIVERGED_INNER                = -7,  /* inner solve failed */
33158c775ebSBarry Smith   SNES_DIVERGED_LOCAL_MIN            = -8,  /* || J^T b || is small, implies converged to local minimum of F() */
332e37c518bSBarry Smith   SNES_DIVERGED_DTOL                 = -9,  /* || F || > divtol*||F_initial|| */
33307b62357SFande Kong   SNES_DIVERGED_JACOBIAN_DOMAIN      = -10, /* Jacobian calculation does not make sense */
3341c6b2ff8SBarry Smith   SNES_DIVERGED_TR_DELTA             = -11,
335c78072a7SPatrick Sanan   SNES_CONVERGED_TR_DELTA_DEPRECATED = -11,
336379ef8d2SAlexander   SNES_DIVERGED_USER                 = -12, /* The user has indicated divergence for an arbitrary reason */
337*76c63389SBarry Smith   SNES_DIVERGED_OBJECTIVE_DOMAIN     = -13,
338*76c63389SBarry Smith   SNES_DIVERGED_OBJECTIVE_NANORINF   = -14,
339e37c518bSBarry Smith 
3409371c9d4SSatish Balay   SNES_CONVERGED_ITERATING = 0
3419371c9d4SSatish Balay } SNESConvergedReason;
342014dd563SJed Brown PETSC_EXTERN const char *const *SNESConvergedReasons;
343184914b5SBarry Smith 
344c838673bSBarry Smith /*MC
345*76c63389SBarry Smith    SNES_CONVERGED_FNORM_ABS - $||F|| \le abstol$
346c838673bSBarry Smith 
347c838673bSBarry Smith    Level: beginner
348c838673bSBarry Smith 
3491cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
350c838673bSBarry Smith M*/
351c838673bSBarry Smith 
352c838673bSBarry Smith /*MC
353*76c63389SBarry Smith    SNES_CONVERGED_FNORM_RELATIVE - $||F|| \le rtol*||F(x_0)||$ where $x_0$ is the initial guess
354c838673bSBarry Smith 
355c838673bSBarry Smith    Level: beginner
356c838673bSBarry Smith 
3571cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
358c838673bSBarry Smith M*/
359c838673bSBarry Smith 
360c838673bSBarry Smith /*MC
361*76c63389SBarry Smith   SNES_CONVERGED_SNORM_RELATIVE - The 2-norm of the last step $\le stol * ||x||$ where $x$ is the current
362*76c63389SBarry Smith   solution and $stol$ is the 4th argument to `SNESSetTolerances()`
363c838673bSBarry Smith 
364af27ebaaSBarry Smith   Options Database Key:
365ae9be289SBarry Smith   -snes_stol <stol> - the step tolerance
366ae9be289SBarry Smith 
367c838673bSBarry Smith    Level: beginner
368c838673bSBarry Smith 
3691cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
370c838673bSBarry Smith M*/
371c838673bSBarry Smith 
372c838673bSBarry Smith /*MC
373c838673bSBarry Smith    SNES_DIVERGED_FUNCTION_COUNT - The user provided function has been called more times then the final
37487497f52SBarry Smith    argument to `SNESSetTolerances()`
375c838673bSBarry Smith 
376c838673bSBarry Smith    Level: beginner
377c838673bSBarry Smith 
3781cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
379c838673bSBarry Smith M*/
380c838673bSBarry Smith 
381c838673bSBarry Smith /*MC
38287497f52SBarry Smith    SNES_DIVERGED_DTOL - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
383e37c518bSBarry Smith 
384e37c518bSBarry Smith    Level: beginner
385e37c518bSBarry Smith 
3861cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESSetDivergenceTolerance()`
387e37c518bSBarry Smith M*/
388e37c518bSBarry Smith 
389e37c518bSBarry Smith /*MC
390*76c63389SBarry Smith    SNES_DIVERGED_FUNCTION_NANORINF - the 2-norm of the current function evaluation is not-a-number (NaN) or infinity, this
391*76c63389SBarry Smith    is usually caused by a division of 0 by 0, or infinity.  See `SNESSetFunctionDomainError()`
392*76c63389SBarry Smith 
393*76c63389SBarry Smith    Level: beginner
394*76c63389SBarry Smith 
395*76c63389SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
396*76c63389SBarry Smith M*/
397*76c63389SBarry Smith 
398*76c63389SBarry Smith /*MC
399*76c63389SBarry Smith    SNES_DIVERGED_FUNCTION_DOMAIN - the function provided with `SNESSetFunction()` called `SNESSetFunctionDomainError()` and
400*76c63389SBarry Smith    the solver could not recoverer.
401*76c63389SBarry Smith 
402*76c63389SBarry Smith    Level: beginner
403*76c63389SBarry Smith 
404*76c63389SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
405*76c63389SBarry Smith M*/
406*76c63389SBarry Smith 
407*76c63389SBarry Smith /*MC
408*76c63389SBarry Smith    SNES_DIVERGED_OBJECTIVE_DOMAIN - the function provided with `SNESSetObjective()` called `SNESSetObjectiveDomainError()` and
409*76c63389SBarry Smith    the solver could not recoverer.
410*76c63389SBarry Smith 
411*76c63389SBarry Smith    Level: beginner
412*76c63389SBarry Smith 
413*76c63389SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
414*76c63389SBarry Smith M*/
415*76c63389SBarry Smith 
416*76c63389SBarry Smith /*MC
417*76c63389SBarry Smith    SNES_DIVERGED_JACOBIAN_DOMAIN - the function provided with `SNESSetJacobian()` called `SNESSetJacobianDomainError()`
418c838673bSBarry Smith 
419c838673bSBarry Smith    Level: beginner
420c838673bSBarry Smith 
4211cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
422c838673bSBarry Smith M*/
423c838673bSBarry Smith 
424c838673bSBarry Smith /*MC
425c838673bSBarry Smith    SNES_DIVERGED_MAX_IT - SNESSolve() has reached the maximum number of iterations requested
426c838673bSBarry Smith 
427c838673bSBarry Smith    Level: beginner
428c838673bSBarry Smith 
4291cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
430c838673bSBarry Smith M*/
431c838673bSBarry Smith 
432c838673bSBarry Smith /*MC
43387497f52SBarry Smith    SNES_DIVERGED_LINE_SEARCH - The line search has failed. This only occurs for a `SNES` solvers that use a line search
434c838673bSBarry Smith 
435c838673bSBarry Smith    Level: beginner
436c838673bSBarry Smith 
4371cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESLineSearch`
438c838673bSBarry Smith M*/
439c838673bSBarry Smith 
440c838673bSBarry Smith /*MC
44146a9e3ceSBarry Smith    SNES_DIVERGED_LOCAL_MIN - the algorithm seems to have stagnated at a local minimum that is not zero.
44287497f52SBarry Smith    See the manual page for `SNESConvergedReason` for more details
443c838673bSBarry Smith 
444c838673bSBarry Smith    Level: beginner
445c838673bSBarry Smith 
4461cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
447c838673bSBarry Smith M*/
448c838673bSBarry Smith 
449c838673bSBarry Smith /*MC
45087497f52SBarry Smith    SNES_CONERGED_ITERATING - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
451c838673bSBarry Smith 
452c838673bSBarry Smith    Level: beginner
453c838673bSBarry Smith 
4541cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
455c838673bSBarry Smith M*/
456c838673bSBarry Smith 
45712651944SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetConvergenceTest(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *), void *, PetscCtxDestroyFn *);
4588d359177SBarry Smith PETSC_EXTERN PetscErrorCode SNESConvergedDefault(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
459e2a6519dSDmitry Karpeev PETSC_EXTERN PetscErrorCode SNESConvergedSkip(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
460649ef022SMatthew Knepley PETSC_EXTERN PetscErrorCode SNESConvergedCorrectPressure(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetConvergedReason(SNES, SNESConvergedReason *);
462c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode SNESGetConvergedReasonString(SNES, const char **);
46333866048SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetConvergedReason(SNES, SNESConvergedReason);
464ddbbdb52SLois Curfman McInnes 
465edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "SNESConvergedSkip()", ) static inline void SNESSkipConverged(void)
466d71ae5a4SJacob Faibussowitsch { /* never called */
4679371c9d4SSatish Balay }
4688ea1b3e6SJed Brown #define SNESSkipConverged (SNESSkipConverged, SNESConvergedSkip)
46911f088b5SMatthew G Knepley 
4709bcc50f1SBarry Smith /*S
4718434afd1SBarry Smith   SNESInitialGuessFn - A prototype of a `SNES` compute initial guess function that would be passed to `SNESSetComputeInitialGuess()`
4729bcc50f1SBarry Smith 
4739bcc50f1SBarry Smith   Calling Sequence:
4749bcc50f1SBarry Smith + snes  - `SNES` context
4759bcc50f1SBarry Smith . u   - output vector to contain initial guess
4769bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4779bcc50f1SBarry Smith 
4789bcc50f1SBarry Smith   Level: beginner
4799bcc50f1SBarry Smith 
4808434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetComputeInitialGuess()`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESFunctionFn`
4819bcc50f1SBarry Smith S*/
48221bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESInitialGuessFn(SNES snes, Vec u, void *ctx);
4839bcc50f1SBarry Smith 
4849bcc50f1SBarry Smith /*S
4858434afd1SBarry Smith   SNESFunctionFn - A prototype of a `SNES` evaluation function that would be passed to `SNESSetFunction()`
4869bcc50f1SBarry Smith 
4879bcc50f1SBarry Smith   Calling Sequence:
4889bcc50f1SBarry Smith + snes  - `SNES` context
4899bcc50f1SBarry Smith . u   - input vector
4909bcc50f1SBarry Smith . F   - function vector
4919bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4929bcc50f1SBarry Smith 
4939bcc50f1SBarry Smith   Level: beginner
4949bcc50f1SBarry Smith 
4958434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
4969bcc50f1SBarry Smith S*/
49721bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESFunctionFn(SNES snes, Vec u, Vec F, void *ctx);
4989bcc50f1SBarry Smith 
4999bcc50f1SBarry Smith /*S
5008434afd1SBarry Smith   SNESObjectiveFn - A prototype of a `SNES` objective evaluation function that would be passed to `SNESSetObjective()`
5019bcc50f1SBarry Smith 
5029bcc50f1SBarry Smith   Calling Sequence:
5039bcc50f1SBarry Smith + snes  - `SNES` context
5049bcc50f1SBarry Smith . u   - input vector
5059bcc50f1SBarry Smith . o   - output value
5069bcc50f1SBarry Smith - ctx - [optional] user-defined function context
5079bcc50f1SBarry Smith 
5089bcc50f1SBarry Smith   Level: beginner
5099bcc50f1SBarry Smith 
5108434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
5119bcc50f1SBarry Smith S*/
51221bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESObjectiveFn(SNES snes, Vec u, PetscReal *o, void *ctx);
5139bcc50f1SBarry Smith 
5149bcc50f1SBarry Smith /*S
5158434afd1SBarry Smith   SNESJacobianFn - A prototype of a `SNES` Jacobian evaluation function that would be passed to `SNESSetJacobian()`
5169bcc50f1SBarry Smith 
5179bcc50f1SBarry Smith   Calling Sequence:
5189bcc50f1SBarry Smith + snes   - the `SNES` context obtained from `SNESCreate()`
5199bcc50f1SBarry Smith . u    - input vector
5209bcc50f1SBarry Smith . Amat - (approximate) Jacobian matrix
5219bcc50f1SBarry Smith . Pmat - matrix used to construct the preconditioner, often the same as `Amat`
5229bcc50f1SBarry Smith - ctx  - [optional] user-defined context for matrix evaluation routine
5239bcc50f1SBarry Smith 
5249bcc50f1SBarry Smith   Level: beginner
5259bcc50f1SBarry Smith 
5268434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESNGSFn`
5279bcc50f1SBarry Smith S*/
52821bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESJacobianFn(SNES snes, Vec u, Mat Amat, Mat Pmat, void *ctx);
5299bcc50f1SBarry Smith 
5309bcc50f1SBarry Smith /*S
5318434afd1SBarry Smith   SNESNGSFn - A prototype of a `SNES` nonlinear Gauss-Seidel function that would be passed to `SNESSetNGS()`
5329bcc50f1SBarry Smith 
5339bcc50f1SBarry Smith   Calling Sequence:
5349bcc50f1SBarry Smith + snes   - the `SNES` context obtained from `SNESCreate()`
5359bcc50f1SBarry Smith . u    - the current solution, updated in place
536dd8e379bSPierre Jolivet . b    - the right-hand side vector (which may be `NULL`)
5379bcc50f1SBarry Smith - ctx  - [optional] user-defined context for matrix evaluation routine
5389bcc50f1SBarry Smith 
5399bcc50f1SBarry Smith   Level: beginner
5409bcc50f1SBarry Smith 
5418434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`
5429bcc50f1SBarry Smith S*/
54321bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESNGSFn(SNES snes, Vec u, Vec b, void *ctx);
5449bcc50f1SBarry Smith 
54553e5d35bSStefano Zampini /*S
54653e5d35bSStefano Zampini   SNESUpdateFn - A prototype of a `SNES` update function that would be passed to `SNESSetUpdate()`
54753e5d35bSStefano Zampini 
54853e5d35bSStefano Zampini   Calling Sequence:
54953e5d35bSStefano Zampini + snes - `SNES` context
55053e5d35bSStefano Zampini - step - the current iteration index
55153e5d35bSStefano Zampini 
55253e5d35bSStefano Zampini   Level: advanced
55353e5d35bSStefano Zampini 
55453e5d35bSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESSetUpdate()`
55553e5d35bSStefano Zampini S*/
55621bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESUpdateFn(SNES snes, PetscInt step);
55753e5d35bSStefano Zampini 
558b67197daSBarry Smith /* --------- Solving systems of nonlinear equations --------------- */
5598434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetFunction(SNES, Vec, SNESFunctionFn *, void *);
5608434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetFunction(SNES, Vec *, SNESFunctionFn **, void **);
561014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESComputeFunction(SNES, Vec, Vec);
562bbc1464cSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeMFFunction(SNES, Vec, Vec);
56325acbd8eSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESSetInitialFunction(SNES, Vec);
56432f3f7c2SPeter Brune 
5658434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetJacobian(SNES, Mat, Mat, SNESJacobianFn *, void *);
5668434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetJacobian(SNES, Mat *, Mat *, SNESJacobianFn **, void **);
5678434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESObjectiveComputeFunctionDefaultFD;
5688434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefault;
5698434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefaultColor;
57078e7fe0eSHong Zhang PETSC_EXTERN PetscErrorCode SNESPruneJacobianColor(SNES, Mat, Mat);
5718434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetComputeInitialGuess(SNES, SNESInitialGuessFn *, void *);
5728434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetPicard(SNES, Vec, SNESFunctionFn *, Mat, Mat, SNESJacobianFn *, void *);
5738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetPicard(SNES, Vec *, SNESFunctionFn **, Mat *, Mat *, SNESJacobianFn **, void **);
5748434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeFunction;
5758434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeMFFunction;
5768434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESPicardComputeJacobian;
57717bae607SBarry Smith 
5788434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetObjective(SNES, SNESObjectiveFn *, void *);
5798434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetObjective(SNES, SNESObjectiveFn **, void **);
5802a4ee8f2SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeObjective(SNES, Vec, PetscReal *);
5812a4ee8f2SPeter Brune 
58253e5d35bSStefano Zampini PETSC_EXTERN PetscErrorCode SNESSetUpdate(SNES, SNESUpdateFn *);
58353e5d35bSStefano Zampini 
584534ebe21SPeter Brune /*E
58516a05f60SBarry Smith    SNESNormSchedule - Frequency with which the norm is computed during a nonliner solve
586534ebe21SPeter Brune 
587dc4c0fb0SBarry Smith    Values:
588dc4c0fb0SBarry Smith +   `SNES_NORM_DEFAULT`            - use the default behavior for the current `SNESType`
589dc4c0fb0SBarry Smith .   `SNES_NORM_NONE`               - avoid all norm computations
590dc4c0fb0SBarry Smith .   `SNES_NORM_ALWAYS`             - compute the norms whenever possible
591dc4c0fb0SBarry Smith .   `SNES_NORM_INITIAL_ONLY`       - compute the norm only when the algorithm starts
592dc4c0fb0SBarry Smith .   `SNES_NORM_FINAL_ONLY`         - compute the norm only when the algorithm finishes
593dc4c0fb0SBarry Smith -   `SNES_NORM_INITIAL_FINAL_ONLY` - compute the norm at the start and end of the algorithm
594534ebe21SPeter Brune 
59516a05f60SBarry Smith    Level: advanced
59616a05f60SBarry Smith 
597534ebe21SPeter Brune    Notes:
598dc4c0fb0SBarry Smith    Support for these is highly dependent on the solver.
599dc4c0fb0SBarry Smith 
600dc4c0fb0SBarry Smith    Some options limit the convergence tests that can be used.
601dc4c0fb0SBarry Smith 
602dc4c0fb0SBarry Smith    The `SNES_NORM_NONE` option is most commonly used when the nonlinear solver is being used as a smoother, for example for `SNESFAS`
603dc4c0fb0SBarry Smith 
604534ebe21SPeter Brune    This is primarily used to turn off extra norm and function computation
605534ebe21SPeter Brune    when the solvers are composed.
606534ebe21SPeter Brune 
6071cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
608db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
609534ebe21SPeter Brune E*/
6109371c9d4SSatish Balay typedef enum {
6119371c9d4SSatish Balay   SNES_NORM_DEFAULT            = -1,
612fdacfa88SPeter Brune   SNES_NORM_NONE               = 0,
613365a6726SPeter Brune   SNES_NORM_ALWAYS             = 1,
614fdacfa88SPeter Brune   SNES_NORM_INITIAL_ONLY       = 2,
615fdacfa88SPeter Brune   SNES_NORM_FINAL_ONLY         = 3,
6169371c9d4SSatish Balay   SNES_NORM_INITIAL_FINAL_ONLY = 4
6179371c9d4SSatish Balay } SNESNormSchedule;
618365a6726SPeter Brune PETSC_EXTERN const char *const *const SNESNormSchedules;
6191957e957SBarry Smith 
620534ebe21SPeter Brune /*MC
62116a05f60SBarry Smith    SNES_NORM_NONE - Don't compute function and its L2 norm when possible
622534ebe21SPeter Brune 
623534ebe21SPeter Brune    Level: advanced
624534ebe21SPeter Brune 
625dc4c0fb0SBarry Smith    Note:
626534ebe21SPeter Brune    This is most useful for stationary solvers with a fixed number of iterations used as smoothers.
627534ebe21SPeter Brune 
6281cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_DEFAULT`
629534ebe21SPeter Brune M*/
630534ebe21SPeter Brune 
631fdacfa88SPeter Brune /*MC
632365a6726SPeter Brune    SNES_NORM_ALWAYS - Compute the function and its L2 norm at each iteration.
633534ebe21SPeter Brune 
634fdacfa88SPeter Brune    Level: advanced
635fdacfa88SPeter Brune 
636dc4c0fb0SBarry Smith    Note:
637fdacfa88SPeter Brune    Most solvers will use this no matter what norm type is passed to them.
638fdacfa88SPeter Brune 
6391cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_NONE`
640fdacfa88SPeter Brune M*/
641534ebe21SPeter Brune 
642534ebe21SPeter Brune /*MC
643534ebe21SPeter Brune    SNES_NORM_INITIAL_ONLY - Compute the function and its L2 at iteration 0, but do not update it.
644534ebe21SPeter Brune 
645534ebe21SPeter Brune    Level: advanced
646534ebe21SPeter Brune 
647534ebe21SPeter Brune    Notes:
648dc4c0fb0SBarry Smith    This method is useful in composed methods, when a true solution might actually be found before `SNESSolve()` is called.
649534ebe21SPeter Brune    This option enables the solve to abort on the zeroth iteration if this is the case.
650534ebe21SPeter Brune 
651534ebe21SPeter Brune    For solvers that require the computation of the L2 norm of the function as part of the method, this merely cancels
652534ebe21SPeter Brune    the norm computation at the last iteration (if possible).
653534ebe21SPeter Brune 
6541cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_FINAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
655534ebe21SPeter Brune M*/
656534ebe21SPeter Brune 
657534ebe21SPeter Brune /*MC
658534ebe21SPeter Brune    SNES_NORM_FINAL_ONLY - Compute the function and its L2 norm on only the final iteration.
659534ebe21SPeter Brune 
660534ebe21SPeter Brune    Level: advanced
661534ebe21SPeter Brune 
662dc4c0fb0SBarry Smith    Note:
663534ebe21SPeter Brune    For solvers that require the computation of the L2 norm of the function as part of the method, behaves
66487497f52SBarry Smith    exactly as `SNES_NORM_DEFAULT`.  This method is useful when the function is gotten after `SNESSolve()` and
665534ebe21SPeter Brune    used in subsequent computation for methods that do not need the norm computed during the rest of the
666534ebe21SPeter Brune    solution procedure.
667534ebe21SPeter Brune 
6681cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_INITIAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
669534ebe21SPeter Brune M*/
670534ebe21SPeter Brune 
671534ebe21SPeter Brune /*MC
672534ebe21SPeter Brune    SNES_NORM_INITIAL_FINAL_ONLY - Compute the function and its L2 norm on only the initial and final iterations.
673534ebe21SPeter Brune 
674534ebe21SPeter Brune    Level: advanced
675534ebe21SPeter Brune 
676dc4c0fb0SBarry Smith    Note:
67787497f52SBarry Smith    This method combines the benefits of `SNES_NORM_INITIAL_ONLY` and `SNES_NORM_FINAL_ONLY`.
678534ebe21SPeter Brune 
6791cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_SNES_NORM_INITIAL_ONLY`, `SNES_NORM_FINAL_ONLY`
680534ebe21SPeter Brune M*/
681534ebe21SPeter Brune 
682365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetNormSchedule(SNES, SNESNormSchedule);
683365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetNormSchedule(SNES, SNESNormSchedule *);
684c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetFunctionNorm(SNES, PetscReal);
685c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESGetFunctionNorm(SNES, PetscReal *);
686c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetUpdateNorm(SNES, PetscReal *);
687c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetSolutionNorm(SNES, PetscReal *);
688534ebe21SPeter Brune 
68947073ea2SPeter Brune /*E
69047073ea2SPeter Brune    SNESFunctionType - Type of function computed
69147073ea2SPeter Brune 
692dc4c0fb0SBarry Smith    Values:
693dc4c0fb0SBarry Smith +  `SNES_FUNCTION_DEFAULT`          - the default behavior for the current `SNESType`
694dc4c0fb0SBarry Smith .  `SNES_FUNCTION_UNPRECONDITIONED` - the original function provided
695dc4c0fb0SBarry Smith -  `SNES_FUNCTION_PRECONDITIONED`   - the modification of the function by the preconditioner
69647073ea2SPeter Brune 
69716a05f60SBarry Smith    Level: advanced
69816a05f60SBarry Smith 
699dc4c0fb0SBarry Smith    Note:
700dc4c0fb0SBarry Smith    Support for these is dependent on the solver.
701dc4c0fb0SBarry Smith 
7021cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
703db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
70447073ea2SPeter Brune E*/
7059371c9d4SSatish Balay typedef enum {
7069371c9d4SSatish Balay   SNES_FUNCTION_DEFAULT          = -1,
70747073ea2SPeter Brune   SNES_FUNCTION_UNPRECONDITIONED = 0,
7089371c9d4SSatish Balay   SNES_FUNCTION_PRECONDITIONED   = 1
7099371c9d4SSatish Balay } SNESFunctionType;
71047073ea2SPeter Brune PETSC_EXTERN const char *const *const SNESFunctionTypes;
71147073ea2SPeter Brune 
71247073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetFunctionType(SNES, SNESFunctionType);
71347073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetFunctionType(SNES, SNESFunctionType *);
714f1c6b773SPeter Brune 
7158434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNGS(SNES, SNESNGSFn *, void *);
7168434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNGS(SNES, SNESNGSFn **, void **);
717be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGS(SNES, Vec, Vec);
718b6266c6eSPeter Brune 
719be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetSweeps(SNES, PetscInt);
720be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetSweeps(SNES, PetscInt *);
721be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt);
722be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
723badc63e7SPeter Brune 
7244fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES, PetscBool);
7254fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES, PetscBool *);
7264fc747eaSLawrence Mitchell 
7273ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode SNESShellGetContext(SNES, void *);
728014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESShellSetContext(SNES, void *);
729014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESShellSetSolve(SNES, PetscErrorCode (*)(SNES, Vec));
730646217ecSPeter Brune 
731c5ae4b9aSBarry Smith /* --------- Routines specifically for line search methods --------------- */
732c5ae4b9aSBarry Smith 
733872b6db9SPeter Brune /*S
73487497f52SBarry Smith    SNESLineSearch - Abstract PETSc object that manages line-search operations for nonlinear solvers
7359e764e56SPeter Brune 
7369e764e56SPeter Brune    Level: beginner
7379e764e56SPeter Brune 
7381cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNES`
7399e764e56SPeter Brune S*/
740907376e6SBarry Smith typedef struct _p_LineSearch *SNESLineSearch;
7419e764e56SPeter Brune 
7429e764e56SPeter Brune /*J
7430b4b7b1cSBarry Smith    SNESLineSearchType - String with the name of a PETSc line search method `SNESLineSearch`. Provides all the linesearches for the nonlinear solvers, `SNES`,
7440b4b7b1cSBarry Smith                         in PETSc.
7459e764e56SPeter Brune 
746ceaaa498SBarry Smith    Values:
747ceaaa498SBarry Smith +  `SNESLINESEARCHBASIC`     - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
748a99ef635SJonas Heinzmann .  `SNESLINESEARCHBT`        - Backtracking line search over the L2 norm of the function or an objective function
749a99ef635SJonas Heinzmann .  `SNESLINESEARCHSECANT`    - Secant line search over the L2 norm of the function or an objective function
7500b4b7b1cSBarry Smith .  `SNESLINESEARCHCP`        - Critical point secant line search assuming $F(x) = \nabla G(x)$ for some unknown $G(x)$
751ceaaa498SBarry Smith .  `SNESLINESEARCHNLEQERR`   - Affine-covariant error-oriented linesearch
752a99ef635SJonas Heinzmann -  `SNESLINESEARCHBISECTION` - bisection line search for a root in the directional derivative
753ceaaa498SBarry Smith -  `SNESLINESEARCHSHELL`     - User provided `SNESLineSearch` implementation
754ceaaa498SBarry Smith 
75595bd0b28SBarry Smith    Level: beginner
75695bd0b28SBarry Smith 
7570b4b7b1cSBarry Smith    Note:
7580b4b7b1cSBarry Smith    Use `SNESLineSearchSetType()` or the options database key `-snes_linesearch_type` to set
7590b4b7b1cSBarry Smith    the specific line search algorithm to use with a given `SNES` object. Not all `SNESType` can utilize a line search.
7600b4b7b1cSBarry Smith 
7611cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetType()`, `SNES`
7629e764e56SPeter Brune J*/
76319fd82e9SBarry Smith typedef const char *SNESLineSearchType;
7649e764e56SPeter Brune #define SNESLINESEARCHBT        "bt"
765d4c6564cSPatrick Farrell #define SNESLINESEARCHNLEQERR   "nleqerr"
766c87759e9SPeter Brune #define SNESLINESEARCHBASIC     "basic"
7670b00b554SBarry Smith #define SNESLINESEARCHNONE      "none"
768a99ef635SJonas Heinzmann #define SNESLINESEARCHSECANT    "secant"
769392273beSMatthew G. Knepley #define SNESLINESEARCHL2        PETSC_DEPRECATED_MACRO(3, 24, 0, "SNESLINESEARCHSECANT", ) "secant"
7709e764e56SPeter Brune #define SNESLINESEARCHCP        "cp"
7719e764e56SPeter Brune #define SNESLINESEARCHSHELL     "shell"
772b5badacbSBarry Smith #define SNESLINESEARCHNCGLINEAR "ncglinear"
773b9d635d7SJonas Heinzmann #define SNESLINESEARCHBISECTION "bisection"
7749e764e56SPeter Brune 
775140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESList;
776014dd563SJed Brown PETSC_EXTERN PetscClassId      SNESLINESEARCH_CLASSID;
777140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESLineSearchList;
7789e764e56SPeter Brune 
779b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_LINEAR    1
780b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_QUADRATIC 2
781b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_CUBIC     3
7829e764e56SPeter Brune 
7839bcc50f1SBarry Smith /*S
7848434afd1SBarry Smith   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that projects a vector onto the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
7859bcc50f1SBarry Smith 
7869bcc50f1SBarry Smith   Calling Sequence:
7879bcc50f1SBarry Smith + snes  - `SNES` context
7889bcc50f1SBarry Smith - u     - the vector to project to the bounds
7899bcc50f1SBarry Smith 
7909bcc50f1SBarry Smith   Level: advanced
7919bcc50f1SBarry Smith 
7929bcc50f1SBarry Smith   Note:
7938434afd1SBarry Smith   The deprecated `SNESLineSearchVIProjectFunc` still works as a replacement for `SNESLineSearchVIProjectFn` *.
7949bcc50f1SBarry Smith 
7959bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
7969bcc50f1SBarry Smith S*/
79721bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode             SNESLineSearchVIProjectFn(SNES snes, Vec u);
79821bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVIProjectFn *SNESLineSearchVIProjectFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVIProjectFn*", );
7999bcc50f1SBarry Smith 
8009bcc50f1SBarry Smith /*S
8018434afd1SBarry Smith   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that computes the norm of the active set variables in a vector in a VI solve,
8029bcc50f1SBarry Smith   passed to `SNESLineSearchSetVIFunctions()`
8039bcc50f1SBarry Smith 
8049bcc50f1SBarry Smith   Calling Sequence:
8059bcc50f1SBarry Smith + snes  - `SNES` context
8069bcc50f1SBarry Smith . f     - the vector to compute the norm of
8079bcc50f1SBarry Smith . u     - the current solution, entries that are on the VI bounds are ignored
8089bcc50f1SBarry Smith - fnorm - the resulting norm
8099bcc50f1SBarry Smith 
8109bcc50f1SBarry Smith   Level: advanced
8119bcc50f1SBarry Smith 
8129bcc50f1SBarry Smith   Note:
8138434afd1SBarry Smith   The deprecated `SNESLineSearchVINormFunc` still works as a replacement for `SNESLineSearchVINormFn` *.
8149bcc50f1SBarry Smith 
8159bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
8169bcc50f1SBarry Smith S*/
81721bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode          SNESLineSearchVINormFn(SNES snes, Vec f, Vec u, PetscReal *fnorm);
81821bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVINormFn *SNESLineSearchVINormFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVINormFnn*", );
8199bcc50f1SBarry Smith 
820d5def619SJonas Heinzmann /*S
821d5def619SJonas Heinzmann   SNESLineSearchVIDirDerivFn - A prototype of a `SNES` function that computes the directional derivative considering the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
822d5def619SJonas Heinzmann 
823d5def619SJonas Heinzmann   Calling Sequence:
824d5def619SJonas Heinzmann + snes  - `SNES` context
825d5def619SJonas Heinzmann . f     - the function vector to compute the directional derivative with
826d5def619SJonas Heinzmann . u     - the current solution, entries that are on the VI bounds are ignored
827d5def619SJonas Heinzmann . y     - the direction to compute the directional derivative
828d5def619SJonas Heinzmann - fty   - the resulting directional derivative
829d5def619SJonas Heinzmann 
830d5def619SJonas Heinzmann   Level: advanced
831d5def619SJonas Heinzmann 
832d5def619SJonas Heinzmann .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchVIProjectFn`, `SNESLineSearchVIProjectFn`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetVIFunctions()`
833d5def619SJonas Heinzmann S*/
83421bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchVIDirDerivFn(SNES snes, Vec f, Vec u, Vec y, PetscScalar *fty);
835d5def619SJonas Heinzmann 
83621bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchApplyFn(SNESLineSearch);
83721bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchApplyFn      *SNESLineSearchApplyFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
83821bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchShellApplyFn(SNESLineSearch, void *);
83921bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchShellApplyFn *SNESLineSearchUserFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
8409e764e56SPeter Brune 
841014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate(MPI_Comm, SNESLineSearch *);
842014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchReset(SNESLineSearch);
843014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchView(SNESLineSearch, PetscViewer);
844014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *);
845a80ff896SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetType(SNESLineSearch, SNESLineSearchType *);
84619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetType(SNESLineSearch, SNESLineSearchType);
847014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch);
848ed07d7d7SPeter Brune PETSC_EXTERN PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch, PetscErrorCode (*)(SNES, Vec, Vec));
849014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetUp(SNESLineSearch);
850014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchApply(SNESLineSearch, Vec, Vec, PetscReal *, Vec);
851014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch, Vec, Vec, PetscBool *);
852014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *);
853fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch, PetscInt);
8549e764e56SPeter Brune 
8559bd66eb0SPeter Brune /* set the functions for precheck and postcheck */
85686d74e61SPeter Brune 
857f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx);
858f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);
8599e764e56SPeter Brune 
860f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx);
861f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);
8629e764e56SPeter Brune 
8639bd66eb0SPeter Brune /* set the functions for VI-specific line search operations */
8649bd66eb0SPeter Brune 
865d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn *, SNESLineSearchVINormFn *, SNESLineSearchVIDirDerivFn *);
866d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn **, SNESLineSearchVINormFn **, SNESLineSearchVIDirDerivFn **);
8679bd66eb0SPeter Brune 
8689e764e56SPeter Brune /* pointers to the associated SNES in order to be able to get the function evaluation out */
869014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch, SNES);
870014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch, SNES *);
8719e764e56SPeter Brune 
8729e764e56SPeter Brune /* set and get the parameters and vectors */
873014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
874014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscInt);
8759e764e56SPeter Brune 
876014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch, Vec, Vec, PetscBool *, void *);
87786d74e61SPeter Brune 
878014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch, PetscReal *);
879014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch, PetscReal);
8809e764e56SPeter Brune 
881014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch, PetscReal *);
882014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch, PetscReal);
8839e764e56SPeter Brune 
884ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch, PetscInt *);
885ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch, PetscInt);
88659405d5eSPeter Brune 
887422a814eSBarry Smith /*E
888dc4c0fb0SBarry Smith     SNESLineSearchReason - indication if the line search has succeeded or failed and why
889422a814eSBarry Smith 
890dc4c0fb0SBarry Smith   Values:
891dc4c0fb0SBarry Smith +  `SNES_LINESEARCH_SUCCEEDED`              - the line search succeeded
892dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_NANORINF`        - a not a number of infinity appeared in the computions
893*76c63389SBarry Smith .  `SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN` - the function was evaluated outside of its domain, see `SNESSetFunctionDomainError()`
894*76c63389SBarry Smith .  `SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN`- the objective function was evaluated outside of its domain, see `SNESSetObjectiveDomainError()`
895*76c63389SBarry Smith .  `SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN` - the Jacobian was evaluated outside of its domain, see `SNESSetJacobianDomainError()`
896dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_REDUCT`          - the linear search failed to get the requested decrease in its norm or objective
897dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_USER`            - used by `SNESLINESEARCHNLEQERR` to indicate the user changed the search direction inappropriately
898dc4c0fb0SBarry Smith -  `SNES_LINESEARCH_FAILED_FUNCTION`        - indicates the maximum number of function evaluations allowed has been surpassed, `SNESConvergedReason` is also
899dc4c0fb0SBarry Smith                                               set to `SNES_DIVERGED_FUNCTION_COUNT`
900422a814eSBarry Smith 
90116a05f60SBarry Smith    Level: intermediate
90216a05f60SBarry Smith 
903dc4c0fb0SBarry Smith    Developer Note:
904*76c63389SBarry Smith    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*76c63389SBarry Smith    search failures.
906422a814eSBarry Smith 
9071cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`,
908dc4c0fb0SBarry Smith           `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
909422a814eSBarry Smith E*/
9109371c9d4SSatish Balay typedef enum {
9119371c9d4SSatish Balay   SNES_LINESEARCH_SUCCEEDED,
912422a814eSBarry Smith   SNES_LINESEARCH_FAILED_NANORINF,
913*76c63389SBarry Smith   SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN,
914*76c63389SBarry Smith   SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN,
915*76c63389SBarry Smith   SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN,
916d5b43468SJose E. Roman   SNES_LINESEARCH_FAILED_REDUCT, /* INSUFFICIENT REDUCTION */
917e9b602ebSSatish Balay   SNES_LINESEARCH_FAILED_USER,
9189371c9d4SSatish Balay   SNES_LINESEARCH_FAILED_FUNCTION
9199371c9d4SSatish Balay } SNESLineSearchReason;
920422a814eSBarry Smith 
921422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetReason(SNESLineSearch, SNESLineSearchReason *);
922422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetReason(SNESLineSearch, SNESLineSearchReason);
9239e764e56SPeter Brune 
924014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch, Vec *, Vec *, Vec *, Vec *, Vec *);
925014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch, Vec, Vec, Vec, Vec, Vec);
9269e764e56SPeter Brune 
927014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *);
928014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch, PetscReal, PetscReal, PetscReal);
929014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch);
930014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch, PetscBool);
9319e764e56SPeter Brune 
932dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitor(SNESLineSearch);
93349abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, void *), void *, PetscCtxDestroyFn *);
934d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch, const char[], const char[], const char[], PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *));
935dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch);
936dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch, PetscViewer);
937dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch, PetscViewer *);
938d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch, PetscViewerAndFormat *);
9399e764e56SPeter Brune 
940ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch, const char[]);
941ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch, const char *[]);
9429e764e56SPeter Brune 
9439e764e56SPeter Brune /* Shell interface functions */
94421bb954dSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellSetApply(SNESLineSearch, SNESLineSearchShellApplyFn *, void *);
9458434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellGetApply(SNESLineSearch, SNESLineSearchShellApplyFn **, void **);
9469bcc50f1SBarry Smith 
94721bb954dSBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellSetApply()", ) static inline PetscErrorCode SNESLineSearchShellSetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn *f, void *ctx)
9489bcc50f1SBarry Smith {
9499bcc50f1SBarry Smith   return SNESLineSearchShellSetApply(ls, f, ctx);
9509bcc50f1SBarry Smith }
9519bcc50f1SBarry Smith 
95221bb954dSBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellGetApply()", ) static inline PetscErrorCode SNESLineSearchShellGetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn **f, void **ctx)
9539bcc50f1SBarry Smith {
9549bcc50f1SBarry Smith   return SNESLineSearchShellGetApply(ls, f, ctx);
9559bcc50f1SBarry Smith }
9569e764e56SPeter Brune 
9572f4102e2SPeter Brune /* BT interface functions */
958014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch, PetscReal);
959014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch, PetscReal *);
9602f4102e2SPeter Brune 
9619e764e56SPeter Brune /*register line search types */
962bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchRegister(const char[], PetscErrorCode (*)(SNESLineSearch));
9639e764e56SPeter Brune 
964720c9a41SShri Abhyankar /* Routines for VI solver */
965014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetVariableBounds(SNES, Vec, Vec);
966cf836535SBlaise Bourdin PETSC_EXTERN PetscErrorCode SNESVIGetVariableBounds(SNES, Vec *, Vec *);
967014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
968014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetInactiveSet(SNES, IS *);
969014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetActiveSetIS(SNES, Vec, Vec, IS *);
970014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES, Vec, Vec, PetscReal *);
971d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFtY(SNES, Vec, Vec, Vec, PetscScalar *);
972014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetRedundancyCheck(SNES, PetscErrorCode (*)(SNES, IS, IS *, void *), void *);
973f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeMeritFunction(Vec, PetscReal *, PetscReal *);
974f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeFunction(SNES, Vec, Vec, void *);
975f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMSetVI(DM, IS);
976f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMDestroyVI(DM);
977f5ea5bd2SShri Abhyankar 
978014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTestLocalMin(SNES);
979da9b6338SBarry Smith 
980eef9c623SLois Curfman McInnes /* Should this routine be private? */
981d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian(SNES, Vec, Mat, Mat);
982cbf8f02cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESTestJacobian(SNES, PetscReal *, PetscReal *);
983494a190aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESTestFunction(SNES);
984ddbbdb52SLois Curfman McInnes 
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetDM(SNES, DM);
986014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetDM(SNES, DM *);
987be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPC(SNES, SNES);
988be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPC(SNES, SNES *);
9893ad1a0b9SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESHasNPC(SNES, PetscBool *);
990be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESApplyNPC(SNES, Vec, Vec, Vec);
991be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCFunction(SNES, Vec, PetscReal *);
992be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeFunctionDefaultNPC(SNES, Vec, Vec);
993be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPCSide(SNES, PCSide);
994be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCSide(SNES, PCSide *);
9957601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLineSearch(SNES, SNESLineSearch);
9967601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLineSearch(SNES, SNESLineSearch *);
9976c699258SBarry Smith 
998edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESGetLineSearch()", ) static inline PetscErrorCode SNESGetSNESLineSearch(SNES snes, SNESLineSearch *ls)
999d71ae5a4SJacob Faibussowitsch {
10009371c9d4SSatish Balay   return SNESGetLineSearch(snes, ls);
10019371c9d4SSatish Balay }
1002edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESSetLineSearch()", ) static inline PetscErrorCode SNESSetSNESLineSearch(SNES snes, SNESLineSearch ls)
1003d71ae5a4SJacob Faibussowitsch {
10049371c9d4SSatish Balay   return SNESSetLineSearch(snes, ls);
10059371c9d4SSatish Balay }
10068b7b3213SJed Brown 
1007014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetUpMatrices(SNES);
10088434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunction(DM, SNESFunctionFn *, void *);
10098434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetFunction(DM, SNESFunctionFn **, void **);
101049abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionContextDestroy(DM, PetscCtxDestroyFn *);
10118434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetMFFunction(DM, SNESFunctionFn *, void *);
10128434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetNGS(DM, SNESNGSFn *, void *);
10138434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetNGS(DM, SNESNGSFn **, void **);
10148434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobian(DM, SNESJacobianFn *, void *);
10158434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetJacobian(DM, SNESJacobianFn **, void **);
101649abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianContextDestroy(DM, PetscCtxDestroyFn *);
10178434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetPicard(DM, SNESFunctionFn *, SNESJacobianFn *, void *);
10188434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetPicard(DM, SNESFunctionFn **, SNESJacobianFn **, void **);
10198434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetObjective(DM, SNESObjectiveFn *, void *);
10208434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetObjective(DM, SNESObjectiveFn **, void **);
10214dd50a75SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMSNES(DM, DM);
10226cab3a1bSJed Brown 
102321bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionFn(DMDALocalInfo *, void *, void *, void *);
102421bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianFn(DMDALocalInfo *, void *, Mat, Mat, void *);
102521bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveFn(DMDALocalInfo *, void *, PetscReal *, void *);
1026c9d099b5SPeter Brune 
102721bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionVecFn(DMDALocalInfo *, Vec, Vec, void *);
102821bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianVecFn(DMDALocalInfo *, Vec, Mat, Mat, void *);
102921bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveVecFn(DMDALocalInfo *, Vec, PetscReal *, void *);
10305505f7afSJunchao Zhang 
10318434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocal(DM, InsertMode, DMDASNESFunctionFn *, void *);
10328434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocal(DM, DMDASNESJacobianFn *, void *);
10338434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocal(DM, DMDASNESObjectiveFn *, void *);
10348434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetPicardLocal(DM, InsertMode, DMDASNESFunctionFn *, DMDASNESJacobianFn, void *);
10356cab3a1bSJed Brown 
10368434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocalVec(DM, InsertMode, DMDASNESFunctionVecFn *, void *);
10378434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocalVec(DM, DMDASNESJacobianVecFn *, void *);
10388434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocalVec(DM, DMDASNESObjectiveVecFn *, void *);
10395505f7afSJunchao Zhang 
1040bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMSNESSetBoundaryLocal(DM, PetscErrorCode (*)(DM, Vec, void *), void *);
10416493148fSStefano Zampini PETSC_EXTERN PetscErrorCode DMSNESSetObjectiveLocal(DM, PetscErrorCode (*)(DM, Vec, PetscReal *, void *), void *);
1042ff35dfedSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionLocal(DM, PetscErrorCode (*)(DM, Vec, Vec, void *), void *);
1043d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianLocal(DM, PetscErrorCode (*)(DM, Vec, Mat, Mat, void *), void *);
104428d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetBoundaryLocal(DM, PetscErrorCode (**)(DM, Vec, void *), void **);
10456493148fSStefano Zampini PETSC_EXTERN PetscErrorCode DMSNESGetObjectiveLocal(DM, PetscErrorCode (**)(DM, Vec, PetscReal *, void *), void **);
104628d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetFunctionLocal(DM, PetscErrorCode (**)(DM, Vec, Vec, void *), void **);
104728d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetJacobianLocal(DM, PetscErrorCode (**)(DM, Vec, Mat, Mat, void *), void **);
1048ff35dfedSBarry Smith 
1049b665c654SMatthew G Knepley /* Routines for Multiblock solver */
1050014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetFields(SNES, const char[], PetscInt, const PetscInt *);
1051014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetIS(SNES, const char[], IS);
1052014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
1053014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);
10544bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMultiblockGetSubSNES(SNES, PetscInt *, SNES *[]);
1055aeed3662SMatthew G Knepley 
105637e1895aSJed Brown /*J
105787497f52SBarry Smith    SNESMSType - String with the name of a PETSc `SNESMS` method.
105837e1895aSJed Brown 
105937e1895aSJed Brown    Level: intermediate
106037e1895aSJed Brown 
10611cc06b55SBarry Smith .seealso: [](ch_snes), `SNESMS`, `SNESMSGetType()`, `SNESMSSetType()`, `SNES`
106237e1895aSJed Brown J*/
106319fd82e9SBarry Smith typedef const char *SNESMSType;
106437e1895aSJed Brown #define SNESMSM62       "m62"
106537e1895aSJed Brown #define SNESMSEULER     "euler"
1066a97cb6bcSJed Brown #define SNESMSJAMESON83 "jameson83"
10673847c725SLisandro Dalcin #define SNESMSVLTP11    "vltp11"
1068a97cb6bcSJed Brown #define SNESMSVLTP21    "vltp21"
1069a97cb6bcSJed Brown #define SNESMSVLTP31    "vltp31"
1070a97cb6bcSJed Brown #define SNESMSVLTP41    "vltp41"
1071a97cb6bcSJed Brown #define SNESMSVLTP51    "vltp51"
1072a97cb6bcSJed Brown #define SNESMSVLTP61    "vltp61"
107337e1895aSJed Brown 
107419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSRegister(SNESMSType, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[]);
10754bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMSRegisterAll(void);
107657715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetType(SNES, SNESMSType *);
107719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSSetType(SNES, SNESMSType);
107857715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetDamping(SNES, PetscReal *);
107957715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSSetDamping(SNES, PetscReal);
1080014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSFinalizePackage(void);
1081607a6623SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSInitializePackage(void);
1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSRegisterDestroy(void);
108337e1895aSJed Brown 
1084ceaaa498SBarry Smith /*MC
1085ceaaa498SBarry Smith    SNESNGMRESRestartType - the restart approach used by `SNESNGMRES`
108613a62661SPeter Brune 
1087ceaaa498SBarry Smith   Values:
1088ceaaa498SBarry Smith +   `SNES_NGMRES_RESTART_NONE`       - never restart
1089ceaaa498SBarry Smith .   `SNES_NGMRES_RESTART_DIFFERENCE` - restart based upon difference criteria
1090ceaaa498SBarry Smith -   `SNES_NGMRES_RESTART_PERIODIC`   - restart after a fixed number of iterations
1091ceaaa498SBarry Smith 
1092ceaaa498SBarry Smith   Options Database Keys:
1093ceaaa498SBarry Smith + -snes_ngmres_restart_type <difference,periodic,none> - set the restart type
1094af27ebaaSBarry Smith - -snes_ngmres_restart <30>                            - sets the number of iterations before restart for periodic
1095ceaaa498SBarry Smith 
109695bd0b28SBarry Smith    Level: intermediate
109795bd0b28SBarry Smith 
1098ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1099ceaaa498SBarry Smith           `SNESNGMRESGetRestartType()`, `SNESNGMRESSelectType`
1100ceaaa498SBarry Smith M*/
110113a62661SPeter Brune typedef enum {
110213a62661SPeter Brune   SNES_NGMRES_RESTART_NONE       = 0,
110313a62661SPeter Brune   SNES_NGMRES_RESTART_PERIODIC   = 1,
11049371c9d4SSatish Balay   SNES_NGMRES_RESTART_DIFFERENCE = 2
11059371c9d4SSatish Balay } SNESNGMRESRestartType;
11066a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];
110713a62661SPeter Brune 
1108ceaaa498SBarry Smith /*MC
1109ceaaa498SBarry Smith    SNESNGMRESSelectType - the approach used by `SNESNGMRES` to determine how the candidate solution and
1110ceaaa498SBarry Smith   combined solution are used to create the next iterate.
1111ceaaa498SBarry Smith 
1112ceaaa498SBarry Smith    Values:
1113ceaaa498SBarry Smith +   `SNES_NGMRES_SELECT_NONE`       - choose the combined solution all the time
1114ceaaa498SBarry Smith .   `SNES_NGMRES_SELECT_DIFFERENCE` - choose based upon the selection criteria
1115ceaaa498SBarry Smith -   `SNES_NGMRES_SELECT_LINESEARCH` - choose based upon line search combination
1116ceaaa498SBarry Smith 
1117ceaaa498SBarry Smith   Options Database Key:
1118ceaaa498SBarry Smith . -snes_ngmres_select_type<difference,none,linesearch> - select type
1119ceaaa498SBarry Smith 
112095bd0b28SBarry Smith    Level: intermediate
112195bd0b28SBarry Smith 
1122ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1123ceaaa498SBarry Smith           `SNESNGMRESGetRestartType()`, `SNESNGMRESRestartType`
1124ceaaa498SBarry Smith M*/
112513a62661SPeter Brune typedef enum {
112613a62661SPeter Brune   SNES_NGMRES_SELECT_NONE       = 0,
112713a62661SPeter Brune   SNES_NGMRES_SELECT_DIFFERENCE = 1,
11289371c9d4SSatish Balay   SNES_NGMRES_SELECT_LINESEARCH = 2
11299371c9d4SSatish Balay } SNESNGMRESSelectType;
11306a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESSelectTypes[];
113113a62661SPeter Brune 
1132014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartType(SNES, SNESNGMRESRestartType);
1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetSelectType(SNES, SNESNGMRESSelectType);
113423b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartFmRise(SNES, PetscBool);
113523b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESGetRestartFmRise(SNES, PetscBool *);
113613a62661SPeter Brune 
1137ceaaa498SBarry Smith /*MC
1138ceaaa498SBarry Smith    SNESNCGType - the conjugate update approach for `SNESNCG`
11390a844d1aSPeter Brune 
1140ceaaa498SBarry Smith    Values:
1141ceaaa498SBarry Smith +   `SNES_NCG_FR`  - Fletcher-Reeves update
1142ceaaa498SBarry Smith .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update, the default and the only one that tolerates generalized search directions
1143ceaaa498SBarry Smith .   `SNES_NCG_HS`  - Hestenes-Steifel update
1144ceaaa498SBarry Smith .   `SNES_NCG_DY`  - Dai-Yuan update
1145ceaaa498SBarry Smith -   `SNES_NCG_CD`  - Conjugate Descent update
1146ceaaa498SBarry Smith 
1147ceaaa498SBarry Smith   Options Database Key:
1148ceaaa498SBarry Smith . -snes_ncg_type<fr,prp,hs,dy,cd> - select type
1149ceaaa498SBarry Smith 
115095bd0b28SBarry Smith    Level: intermediate
115195bd0b28SBarry Smith 
1152ceaaa498SBarry Smith .seealso: `SNES, `SNESNCG`, `SNESNCGSetType()`
1153ceaaa498SBarry Smith M*/
11540a844d1aSPeter Brune typedef enum {
11550de8b71eSPeter Brune   SNES_NCG_FR  = 0,
11560de8b71eSPeter Brune   SNES_NCG_PRP = 1,
11570de8b71eSPeter Brune   SNES_NCG_HS  = 2,
11580de8b71eSPeter Brune   SNES_NCG_DY  = 3,
11599371c9d4SSatish Balay   SNES_NCG_CD  = 4
11609371c9d4SSatish Balay } SNESNCGType;
11616a6fc655SJed Brown PETSC_EXTERN const char *const SNESNCGTypes[];
11620a844d1aSPeter Brune 
1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNCGSetType(SNES, SNESNCGType);
116413a62661SPeter Brune 
1165ceaaa498SBarry Smith /*MC
1166ceaaa498SBarry Smith    SNESQNScaleType - the scaling type used by `SNESQN`
1167ceaaa498SBarry Smith 
1168ceaaa498SBarry Smith    Values:
1169ceaaa498SBarry Smith +   `SNES_QN_SCALE_NONE`     - don't scale the problem
1170ceaaa498SBarry Smith .   `SNES_QN_SCALE_SCALAR`   - use Shanno scaling
1171ceaaa498SBarry Smith .   `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
1172ceaaa498SBarry Smith -   `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with `SNESSetJacobian()`
1173ceaaa498SBarry Smith                                computed at the first iteration of `SNESQN` and at ever restart.
1174ceaaa498SBarry Smith 
1175ceaaa498SBarry Smith     Options Database Key:
1176ceaaa498SBarry Smith . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
1177ceaaa498SBarry Smith 
117895bd0b28SBarry Smith    Level: intermediate
117995bd0b28SBarry Smith 
1180ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNRestartType`
1181ceaaa498SBarry Smith M*/
11829371c9d4SSatish Balay typedef enum {
11839371c9d4SSatish Balay   SNES_QN_SCALE_DEFAULT  = 0,
11841efc8c45SPeter Brune   SNES_QN_SCALE_NONE     = 1,
118592f76d53SAlp Dener   SNES_QN_SCALE_SCALAR   = 2,
118692f76d53SAlp Dener   SNES_QN_SCALE_DIAGONAL = 3,
11879371c9d4SSatish Balay   SNES_QN_SCALE_JACOBIAN = 4
11889371c9d4SSatish Balay } SNESQNScaleType;
11896a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNScaleTypes[];
1190ceaaa498SBarry Smith 
1191ceaaa498SBarry Smith /*MC
1192ceaaa498SBarry Smith    SNESQNRestartType - the restart approached used by `SNESQN`
1193ceaaa498SBarry Smith 
1194ceaaa498SBarry Smith    Values:
1195ceaaa498SBarry Smith +   `SNES_QN_RESTART_NONE`     - never restart
1196ceaaa498SBarry Smith .   `SNES_QN_RESTART_POWELL`   - restart based upon descent criteria
1197ceaaa498SBarry Smith -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations
1198ceaaa498SBarry Smith 
1199ceaaa498SBarry Smith   Options Database Keys:
1200ceaaa498SBarry Smith + -snes_qn_restart_type <powell,periodic,none> - set the restart type
1201ceaaa498SBarry Smith - -snes_qn_m <m>                               - sets the number of stored updates and the restart period for periodic
1202ceaaa498SBarry Smith 
120395bd0b28SBarry Smith    Level: intermediate
120495bd0b28SBarry Smith 
1205ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNScaleType`
1206ceaaa498SBarry Smith M*/
12079371c9d4SSatish Balay typedef enum {
12089371c9d4SSatish Balay   SNES_QN_RESTART_DEFAULT  = 0,
12091efc8c45SPeter Brune   SNES_QN_RESTART_NONE     = 1,
12101efc8c45SPeter Brune   SNES_QN_RESTART_POWELL   = 2,
12119371c9d4SSatish Balay   SNES_QN_RESTART_PERIODIC = 3
12129371c9d4SSatish Balay } SNESQNRestartType;
12136a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNRestartTypes[];
1214ceaaa498SBarry Smith 
1215ceaaa498SBarry Smith /*MC
1216ceaaa498SBarry Smith    SNESQNType - the type used by `SNESQN`
1217ceaaa498SBarry Smith 
1218ceaaa498SBarry Smith   Values:
1219ceaaa498SBarry Smith +   `SNES_QN_LBFGS`      - LBFGS variant
1220ceaaa498SBarry Smith .   `SNES_QN_BROYDEN`    - Broyden variant
1221ceaaa498SBarry Smith -   `SNES_QN_BADBROYDEN` - Bad Broyden variant
1222ceaaa498SBarry Smith 
1223ceaaa498SBarry Smith   Options Database Key:
1224ceaaa498SBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
1225ceaaa498SBarry Smith 
122695bd0b28SBarry Smith    Level: intermediate
122795bd0b28SBarry Smith 
1228ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNSetType()`, `SNESQNScaleType`, `SNESQNRestartType`, `SNESQNSetRestartType()`
1229ceaaa498SBarry Smith M*/
12309371c9d4SSatish Balay typedef enum {
12319371c9d4SSatish Balay   SNES_QN_LBFGS      = 0,
12321efc8c45SPeter Brune   SNES_QN_BROYDEN    = 1,
12331efc8c45SPeter Brune   SNES_QN_BADBROYDEN = 2
12341efc8c45SPeter Brune } SNESQNType;
12351efc8c45SPeter Brune PETSC_EXTERN const char *const SNESQNTypes[];
12360c777b0cSPeter Brune 
12371efc8c45SPeter Brune PETSC_EXTERN PetscErrorCode SNESQNSetType(SNES, SNESQNType);
1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetScaleType(SNES, SNESQNScaleType);
1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetRestartType(SNES, SNESQNRestartType);
12400c777b0cSPeter Brune 
1241e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetType(SNES, PCASMType *);
1242e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetType(SNES, PCASMType);
1243ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomains(SNES, PetscInt *, SNES *[], VecScatter *[], VecScatter *[], VecScatter *[]);
1244ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMSetSubdomains(SNES, PetscInt, SNES[], VecScatter[], VecScatter[], VecScatter[]);
1245610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetDamping(SNES, PetscReal);
1246610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetDamping(SNES, PetscReal *);
1247ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomainVecs(SNES, PetscInt *, Vec *[], Vec *[], Vec *[], Vec *[]);
1248d728fb7dSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES, PetscBool);
1249d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetSNES(SNES, PetscInt, SNES *);
1250d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetNumber(SNES, PetscInt *);
1251f10b3e88SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMSetWeight(SNES, Vec);
12520c777b0cSPeter Brune 
1253ad05f2b7SBarry Smith /*E
1254ad05f2b7SBarry Smith   SNESCompositeType - Determines how two or more preconditioners are composed with the `SNESType` of `SNESCOMPOSITE`
1255ad05f2b7SBarry Smith 
1256ad05f2b7SBarry Smith   Values:
1257ad05f2b7SBarry Smith + `SNES_COMPOSITE_ADDITIVE`        - results from application of all preconditioners are added together
1258ad05f2b7SBarry Smith . `SNES_COMPOSITE_MULTIPLICATIVE`  - preconditioners are applied sequentially to the residual freshly
1259ad05f2b7SBarry Smith                                      computed after the previous preconditioner application
1260ad05f2b7SBarry Smith - `SNES_COMPOSITE_ADDITIVEOPTIMAL` - uses a linear combination of the solutions obtained with each preconditioner that approximately minimize the function
1261ad05f2b7SBarry Smith                                      value at the new iteration.
1262ad05f2b7SBarry Smith 
1263ad05f2b7SBarry Smith    Level: beginner
1264ad05f2b7SBarry Smith 
1265ad05f2b7SBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `PCCompositeType`
1266ad05f2b7SBarry Smith E*/
12679371c9d4SSatish Balay typedef enum {
12689371c9d4SSatish Balay   SNES_COMPOSITE_ADDITIVE,
12699371c9d4SSatish Balay   SNES_COMPOSITE_MULTIPLICATIVE,
12709371c9d4SSatish Balay   SNES_COMPOSITE_ADDITIVEOPTIMAL
12719371c9d4SSatish Balay } SNESCompositeType;
1272eed5f15bSPeter Brune PETSC_EXTERN const char *const SNESCompositeTypes[];
1273eed5f15bSPeter Brune 
1274eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetType(SNES, SNESCompositeType);
1275eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeAddSNES(SNES, SNESType);
1276eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeGetSNES(SNES, PetscInt, SNES *);
1277a6b47ab3SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESCompositeGetNumber(SNES, PetscInt *);
12788f626970SPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetDamping(SNES, PetscInt, PetscReal);
1279eed5f15bSPeter Brune 
128068e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetDiscretisationInfo(SNES, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
12814d04e9f1SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetComputeOperator(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
128239fd2e8aSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetComputeFunction(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
128368e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetConstructType(SNES, PCPatchConstructType, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *);
128468e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetCellNumbering(SNES, PetscSection);
1285561742edSMatthew G. Knepley 
1286a055c8caSBarry Smith /*E
1287a055c8caSBarry Smith     SNESFASType - Determines the type of nonlinear multigrid method that is run.
1288a055c8caSBarry Smith 
1289a055c8caSBarry Smith    Values:
129087497f52SBarry Smith +  `SNES_FAS_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `SNESFASSetCycles()`
129187497f52SBarry Smith .  `SNES_FAS_ADDITIVE`                 - additive FAS cycle
129287497f52SBarry Smith .  `SNES_FAS_FULL`                     - full FAS cycle
129387497f52SBarry Smith -  `SNES_FAS_KASKADE`                  - Kaskade FAS cycle
1294a055c8caSBarry Smith 
129516a05f60SBarry Smith    Level: beginner
129616a05f60SBarry Smith 
12971cc06b55SBarry Smith .seealso: [](ch_snes), `SNESFAS`, `PCMGSetType()`, `PCMGType`
1298a055c8caSBarry Smith E*/
12999371c9d4SSatish Balay typedef enum {
13009371c9d4SSatish Balay   SNES_FAS_MULTIPLICATIVE,
13019371c9d4SSatish Balay   SNES_FAS_ADDITIVE,
13029371c9d4SSatish Balay   SNES_FAS_FULL,
13039371c9d4SSatish Balay   SNES_FAS_KASKADE
13049371c9d4SSatish Balay } SNESFASType;
1305a055c8caSBarry Smith PETSC_EXTERN const char *const SNESFASTypes[];
1306a055c8caSBarry Smith 
1307a055c8caSBarry Smith /* called on the finest level FAS instance*/
1308a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetType(SNES, SNESFASType);
1309a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetType(SNES, SNESFASType *);
1310a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLevels(SNES, PetscInt, MPI_Comm *);
1311a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetLevels(SNES, PetscInt *);
1312a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCycleSNES(SNES, PetscInt, SNES *);
1313a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothUp(SNES, PetscInt);
1314a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothDown(SNES, PetscInt);
1315a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetCycles(SNES, PetscInt);
1316d142ab34SLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESFASSetMonitor(SNES, PetscViewerAndFormat *, PetscBool);
1317a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLog(SNES, PetscBool);
1318a055c8caSBarry Smith 
1319a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetGalerkin(SNES, PetscBool);
1320a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetGalerkin(SNES, PetscBool *);
132125acbd8eSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESFASGalerkinFunctionDefault(SNES, Vec, Vec, void *);
1322a055c8caSBarry Smith 
1323a055c8caSBarry Smith /* called on any level -- "Cycle" FAS instance */
1324a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmoother(SNES, SNES *);
1325a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherUp(SNES, SNES *);
1326a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherDown(SNES, SNES *);
1327a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetCorrection(SNES, SNES *);
1328a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInterpolation(SNES, Mat *);
1329a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRestriction(SNES, Mat *);
1330a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInjection(SNES, Mat *);
1331a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRScale(SNES, Vec *);
1332a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleSetCycles(SNES, PetscInt);
1333a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleIsFine(SNES, PetscBool *);
1334a055c8caSBarry Smith 
1335a055c8caSBarry Smith /* called on the (outer) finest level FAS to set/get parameters on any level instance */
1336a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInterpolation(SNES, PetscInt, Mat);
1337a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInterpolation(SNES, PetscInt, Mat *);
1338a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRestriction(SNES, PetscInt, Mat);
1339a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRestriction(SNES, PetscInt, Mat *);
1340a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInjection(SNES, PetscInt, Mat);
1341a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInjection(SNES, PetscInt, Mat *);
1342a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRScale(SNES, PetscInt, Vec);
1343a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRScale(SNES, PetscInt, Vec *);
1344020631bcSPeter Brune PETSC_EXTERN PetscErrorCode SNESFASSetContinuation(SNES, PetscBool);
1345a055c8caSBarry Smith 
1346a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmoother(SNES, PetscInt, SNES *);
1347a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherUp(SNES, PetscInt, SNES *);
1348a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherDown(SNES, PetscInt, SNES *);
1349a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCoarseSolve(SNES, SNES *);
1350a055c8caSBarry Smith 
1351a055c8caSBarry Smith /* parameters for full FAS */
1352a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASFullSetDownSweep(SNES, PetscBool);
1353a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCreateCoarseVec(SNES, Vec *);
1354a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASRestrict(SNES, Vec, Vec);
13556dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullSetTotal(SNES, PetscBool);
13566dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullGetTotal(SNES, PetscBool *);
1357a055c8caSBarry Smith 
13582eace19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetSNESVariableBounds(DM, SNES);
1359f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckDiscretization(SNES, DM, PetscReal, Vec, PetscReal, PetscReal[]);
1360f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckResidual(SNES, DM, Vec, PetscReal, PetscReal *);
1361f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckJacobian(SNES, DM, Vec, PetscReal, PetscBool *, PetscReal *);
13627f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckFromOptions(SNES, Vec);
13638e3a2eefSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESComputeJacobianAction(DM, Vec, Vec, Vec, void *);
13648e3a2eefSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCreateJacobianMF(DM, Vec, void *, Mat *);
136597276fddSZach Atkins 
136697276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALSetFunction(SNES, SNESFunctionFn *, void *ctx);
136797276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALGetFunction(SNES, SNESFunctionFn **, void **ctx);
136897276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALComputeFunction(SNES, Vec, Vec);
136997276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALGetLoadParameter(SNES, PetscReal *);
137097276fddSZach Atkins 
137197276fddSZach Atkins /*MC
137297276fddSZach Atkins    SNESNewtonALCorrectionType - the approach used by `SNESNEWTONAL` to determine
137397276fddSZach Atkins    the correction to the current increment. While the exact correction satisfies
137497276fddSZach Atkins    the constraint surface at every iteration, it also requires solving a quadratic
137597276fddSZach Atkins    equation which may not have real roots. Conversely, the normal correction is more
137697276fddSZach Atkins    efficient and always yields a real correction and is the default.
137797276fddSZach Atkins 
137897276fddSZach Atkins    Values:
137997276fddSZach Atkins +   `SNES_NEWTONAL_CORRECTION_EXACT` - choose the correction which exactly satisfies the constraint
1380d7c1f440SPierre Jolivet -   `SNES_NEWTONAL_CORRECTION_NORMAL` - choose the correction in the updated normal hyper-surface to the constraint surface
138197276fddSZach Atkins 
138297276fddSZach Atkins    Options Database Key:
138397276fddSZach Atkins . -snes_newtonal_correction_type <exact> - select type from <exact,normal>
138497276fddSZach Atkins 
138597276fddSZach Atkins    Level: intermediate
138697276fddSZach Atkins 
138797276fddSZach Atkins .seealso: `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetCorrectionType()`
138897276fddSZach Atkins M*/
138997276fddSZach Atkins typedef enum {
139097276fddSZach Atkins   SNES_NEWTONAL_CORRECTION_EXACT  = 0,
139197276fddSZach Atkins   SNES_NEWTONAL_CORRECTION_NORMAL = 1,
139297276fddSZach Atkins } SNESNewtonALCorrectionType;
139397276fddSZach Atkins PETSC_EXTERN const char *const SNESNewtonALCorrectionTypes[];
139497276fddSZach Atkins 
139597276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALSetCorrectionType(SNES, SNESNewtonALCorrectionType);
1396