xref: /petsc/include/petscsnes.h (revision 21bb954d41cf3cfdec1f910396b278bc0ac10d56)
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);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetFunctionDomainError(SNES, PetscBool *);
249cc6b0f04SFande Kong PETSC_EXTERN PetscErrorCode SNESGetJacobianDomainError(SNES, PetscBool *);
250cc6b0f04SFande Kong PETSC_EXTERN PetscErrorCode SNESSetJacobianDomainError(SNES);
251b351a90bSFande Kong PETSC_EXTERN PetscErrorCode SNESSetCheckJacobianDomainError(SNES, PetscBool);
2528383d7d7SFande Kong PETSC_EXTERN PetscErrorCode SNESGetCheckJacobianDomainError(SNES, PetscBool *);
2536a388c36SPeter Brune 
254edd03b47SJacob Faibussowitsch #define SNES_CONVERGED_TR_DELTA_DEPRECATED SNES_CONVERGED_TR_DELTA PETSC_DEPRECATED_ENUM(3, 12, 0, "SNES_DIVERGED_TR_DELTA", )
255435da068SBarry Smith /*E
256dc4c0fb0SBarry Smith     SNESConvergedReason - reason a `SNESSolve()` was determined to have converged or diverged
257435da068SBarry Smith 
258dc4c0fb0SBarry Smith    Values:
259dc4c0fb0SBarry Smith +  `SNES_CONVERGED_FNORM_ABS`      - 2-norm(F) <= abstol
260dc4c0fb0SBarry Smith .  `SNES_CONVERGED_FNORM_RELATIVE` - 2-norm(F) <= rtol*2-norm(F(x_0)) where x_0 is the initial guess
261dc4c0fb0SBarry Smith .  `SNES_CONVERGED_SNORM_RELATIVE` - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
262379ef8d2SAlexander .  `SNES_CONVERGED_USER`           - The user has indicated convergence for an arbitrary reason
263dc4c0fb0SBarry Smith .  `SNES_DIVERGED_FUNCTION_COUNT`  - The user provided function has been called more times than the maximum set in `SNESSetTolerances()`
264dc4c0fb0SBarry Smith .  `SNES_DIVERGED_DTOL`            - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
265dc4c0fb0SBarry Smith .  `SNES_DIVERGED_FNORM_NAN`       - the 2-norm of the current function evaluation is not-a-number (NaN), this
266dc4c0fb0SBarry Smith                                      is usually caused by a division of 0 by 0.
267dc4c0fb0SBarry Smith .  `SNES_DIVERGED_MAX_IT`          - `SNESSolve()` has reached the maximum number of iterations requested
268dc4c0fb0SBarry Smith .  `SNES_DIVERGED_LINE_SEARCH`     - The line search has failed. This only occurs for `SNES` solvers that use a line search
269dc4c0fb0SBarry Smith .  `SNES_DIVERGED_LOCAL_MIN`       - the algorithm seems to have stagnated at a local minimum that is not zero.
270379ef8d2SAlexander .  `SNES_DIVERGED_USER`            - The user has indicated divergence for an arbitrary reason
271379ef8d2SAlexander -  `SNES_CONVERGED_ITERATING       - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
272f203c74bSBarry Smith 
27316a05f60SBarry Smith    Level: beginner
27416a05f60SBarry Smith 
275dc4c0fb0SBarry Smith     Notes:
2760b4b7b1cSBarry 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
2770b4b7b1cSBarry Smith    (in this case we recommend
278dc4c0fb0SBarry Smith    testing with `-pc_type lu` to eliminate the linear solver as the cause of the problem).
279dc4c0fb0SBarry Smith 
2800b4b7b1cSBarry Smith    `SNES_DIVERGED_LOCAL_MIN` can only occur when using a `SNES` solver that uses a line search (`SNESLineSearch`).
281f203c74bSBarry Smith    The line search wants to minimize Q(alpha) = 1/2 || F(x + alpha s) ||^2_2  this occurs
282f203c74bSBarry 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
283f203c74bSBarry Smith    you get Q'(alpha) = -F(x)^T F'(x)^(-1)^T F'(x+alpha s)F(x+alpha s); when alpha = 0
284f203c74bSBarry Smith    Q'(0) = - ||F(x)||^2_2 which is always NEGATIVE if F'(x) is invertible. This means the Newton
285f203c74bSBarry Smith    direction is a descent direction and the line search should succeed if alpha is small enough.
286f203c74bSBarry Smith 
287f203c74bSBarry Smith    If F'(x) is NOT invertible AND F'(x)^T F(x) = 0 then Q'(0) = 0 and the Newton direction
288f203c74bSBarry Smith    is NOT a descent direction so the line search will fail. All one can do at this point
289f203c74bSBarry Smith    is change the initial guess and try again.
290f203c74bSBarry Smith 
291f203c74bSBarry Smith    An alternative explanation: Newton's method can be regarded as replacing the function with
292f203c74bSBarry Smith    its linear approximation and minimizing the 2-norm of that. That is F(x+s) approx F(x) + F'(x)s
293f203c74bSBarry Smith    so we minimize || F(x) + F'(x) s ||^2_2; do this using Least Squares. If F'(x) is invertible then
294f203c74bSBarry 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
2956aad120cSJose E. Roman    exists a nontrivial (that is F'(x)s != 0) solution to the equation and this direction is
296f203c74bSBarry 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)
297f203c74bSBarry 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)) != 0
298f203c74bSBarry Smith    and F'(x)^T F'(x) has no negative eigenvalues Q'(0) < 0 so s is a descent direction and the line
299f203c74bSBarry Smith    search should succeed for small enough alpha.
300f203c74bSBarry Smith 
301f203c74bSBarry Smith    Note that this RARELY happens in practice. Far more likely the linear system is not being solved
302f203c74bSBarry Smith    (well enough?) or the Jacobian is wrong.
303f203c74bSBarry Smith 
30487497f52SBarry Smith    `SNES_DIVERGED_MAX_IT` means that the solver reached the maximum number of iterations without satisfying any
30587497f52SBarry Smith    convergence criteria. `SNES_CONVERGED_ITS` means that `SNESConvergedSkip()` was chosen as the convergence test;
3064d0a8057SBarry Smith    thus the usual convergence criteria have not been checked and may or may not be satisfied.
307f203c74bSBarry Smith 
3081cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`, `SNESSetTolerances()`
309435da068SBarry Smith E*/
310184914b5SBarry Smith typedef enum {                       /* converged */
31101b82886SBarry Smith   SNES_CONVERGED_FNORM_ABS      = 2, /* ||F|| < atol */
31201b82886SBarry Smith   SNES_CONVERGED_FNORM_RELATIVE = 3, /* ||F|| < rtol*||F_initial|| */
3135358d0d4SBarry Smith   SNES_CONVERGED_SNORM_RELATIVE = 4, /* Newton computed step size small; || delta x || < stol || x || */
3143f149594SLisandro Dalcin   SNES_CONVERGED_ITS            = 5, /* maximum iterations reached */
315ce78bad3SBarry 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 */
316379ef8d2SAlexander   SNES_CONVERGED_USER           = 7, /* The user has indicated convergence for an arbitrary reason */
317184914b5SBarry Smith   /* diverged */
31846a9e3ceSBarry Smith   SNES_DIVERGED_FUNCTION_DOMAIN      = -1, /* the new x location passed the function is not in the domain of F */
319184914b5SBarry Smith   SNES_DIVERGED_FUNCTION_COUNT       = -2,
32046a9e3ceSBarry Smith   SNES_DIVERGED_LINEAR_SOLVE         = -3, /* the linear solve failed */
321184914b5SBarry Smith   SNES_DIVERGED_FNORM_NAN            = -4,
322184914b5SBarry Smith   SNES_DIVERGED_MAX_IT               = -5,
323647a2e1fSBarry Smith   SNES_DIVERGED_LINE_SEARCH          = -6,  /* the line search failed */
3241e633543SBarry Smith   SNES_DIVERGED_INNER                = -7,  /* inner solve failed */
32558c775ebSBarry Smith   SNES_DIVERGED_LOCAL_MIN            = -8,  /* || J^T b || is small, implies converged to local minimum of F() */
326e37c518bSBarry Smith   SNES_DIVERGED_DTOL                 = -9,  /* || F || > divtol*||F_initial|| */
32707b62357SFande Kong   SNES_DIVERGED_JACOBIAN_DOMAIN      = -10, /* Jacobian calculation does not make sense */
3281c6b2ff8SBarry Smith   SNES_DIVERGED_TR_DELTA             = -11,
329c78072a7SPatrick Sanan   SNES_CONVERGED_TR_DELTA_DEPRECATED = -11,
330379ef8d2SAlexander   SNES_DIVERGED_USER                 = -12, /* The user has indicated divergence for an arbitrary reason */
331e37c518bSBarry Smith 
3329371c9d4SSatish Balay   SNES_CONVERGED_ITERATING = 0
3339371c9d4SSatish Balay } SNESConvergedReason;
334014dd563SJed Brown PETSC_EXTERN const char *const *SNESConvergedReasons;
335184914b5SBarry Smith 
336c838673bSBarry Smith /*MC
337c838673bSBarry Smith    SNES_CONVERGED_FNORM_ABS - 2-norm(F) <= abstol
338c838673bSBarry Smith 
339c838673bSBarry Smith    Level: beginner
340c838673bSBarry Smith 
3411cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
342c838673bSBarry Smith M*/
343c838673bSBarry Smith 
344c838673bSBarry Smith /*MC
345c838673bSBarry Smith    SNES_CONVERGED_FNORM_RELATIVE - 2-norm(F) <= rtol*2-norm(F(x_0)) where x_0 is the initial guess
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
353c60f73f4SPeter Brune   SNES_CONVERGED_SNORM_RELATIVE - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
35487497f52SBarry Smith   solution and stol is the 4th argument to `SNESSetTolerances()`
355c838673bSBarry Smith 
356af27ebaaSBarry Smith   Options Database Key:
357ae9be289SBarry Smith   -snes_stol <stol> - the step tolerance
358ae9be289SBarry Smith 
359c838673bSBarry Smith    Level: beginner
360c838673bSBarry Smith 
3611cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
362c838673bSBarry Smith M*/
363c838673bSBarry Smith 
364c838673bSBarry Smith /*MC
365c838673bSBarry Smith    SNES_DIVERGED_FUNCTION_COUNT - The user provided function has been called more times then the final
36687497f52SBarry Smith    argument to `SNESSetTolerances()`
367c838673bSBarry Smith 
368c838673bSBarry Smith    Level: beginner
369c838673bSBarry Smith 
3701cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
371c838673bSBarry Smith M*/
372c838673bSBarry Smith 
373c838673bSBarry Smith /*MC
37487497f52SBarry Smith    SNES_DIVERGED_DTOL - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
375e37c518bSBarry Smith 
376e37c518bSBarry Smith    Level: beginner
377e37c518bSBarry Smith 
3781cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESSetDivergenceTolerance()`
379e37c518bSBarry Smith M*/
380e37c518bSBarry Smith 
381e37c518bSBarry Smith /*MC
382c838673bSBarry Smith    SNES_DIVERGED_FNORM_NAN - the 2-norm of the current function evaluation is not-a-number (NaN), this
383c838673bSBarry Smith    is usually caused by a division of 0 by 0.
384c838673bSBarry Smith 
385c838673bSBarry Smith    Level: beginner
386c838673bSBarry Smith 
3871cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
388c838673bSBarry Smith M*/
389c838673bSBarry Smith 
390c838673bSBarry Smith /*MC
391c838673bSBarry Smith    SNES_DIVERGED_MAX_IT - SNESSolve() has reached the maximum number of iterations requested
392c838673bSBarry Smith 
393c838673bSBarry Smith    Level: beginner
394c838673bSBarry Smith 
3951cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
396c838673bSBarry Smith M*/
397c838673bSBarry Smith 
398c838673bSBarry Smith /*MC
39987497f52SBarry Smith    SNES_DIVERGED_LINE_SEARCH - The line search has failed. This only occurs for a `SNES` solvers that use a line search
400c838673bSBarry Smith 
401c838673bSBarry Smith    Level: beginner
402c838673bSBarry Smith 
4031cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESLineSearch`
404c838673bSBarry Smith M*/
405c838673bSBarry Smith 
406c838673bSBarry Smith /*MC
40746a9e3ceSBarry Smith    SNES_DIVERGED_LOCAL_MIN - the algorithm seems to have stagnated at a local minimum that is not zero.
40887497f52SBarry Smith    See the manual page for `SNESConvergedReason` for more details
409c838673bSBarry Smith 
410c838673bSBarry Smith    Level: beginner
411c838673bSBarry Smith 
4121cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
413c838673bSBarry Smith M*/
414c838673bSBarry Smith 
415c838673bSBarry Smith /*MC
41687497f52SBarry Smith    SNES_CONERGED_ITERATING - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
417c838673bSBarry Smith 
418c838673bSBarry Smith    Level: beginner
419c838673bSBarry Smith 
4201cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
421c838673bSBarry Smith M*/
422c838673bSBarry Smith 
423014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetConvergenceTest(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
4248d359177SBarry Smith PETSC_EXTERN PetscErrorCode SNESConvergedDefault(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
425e2a6519dSDmitry Karpeev PETSC_EXTERN PetscErrorCode SNESConvergedSkip(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
426649ef022SMatthew Knepley PETSC_EXTERN PetscErrorCode SNESConvergedCorrectPressure(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
427014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetConvergedReason(SNES, SNESConvergedReason *);
428c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode SNESGetConvergedReasonString(SNES, const char **);
42933866048SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetConvergedReason(SNES, SNESConvergedReason);
430ddbbdb52SLois Curfman McInnes 
431edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "SNESConvergedSkip()", ) static inline void SNESSkipConverged(void)
432d71ae5a4SJacob Faibussowitsch { /* never called */
4339371c9d4SSatish Balay }
4348ea1b3e6SJed Brown #define SNESSkipConverged (SNESSkipConverged, SNESConvergedSkip)
43511f088b5SMatthew G Knepley 
4369bcc50f1SBarry Smith /*S
4378434afd1SBarry Smith   SNESInitialGuessFn - A prototype of a `SNES` compute initial guess function that would be passed to `SNESSetComputeInitialGuess()`
4389bcc50f1SBarry Smith 
4399bcc50f1SBarry Smith   Calling Sequence:
4409bcc50f1SBarry Smith + snes  - `SNES` context
4419bcc50f1SBarry Smith . u   - output vector to contain initial guess
4429bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4439bcc50f1SBarry Smith 
4449bcc50f1SBarry Smith   Level: beginner
4459bcc50f1SBarry Smith 
4468434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetComputeInitialGuess()`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESFunctionFn`
4479bcc50f1SBarry Smith S*/
448*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESInitialGuessFn(SNES snes, Vec u, void *ctx);
4499bcc50f1SBarry Smith 
4509bcc50f1SBarry Smith /*S
4518434afd1SBarry Smith   SNESFunctionFn - A prototype of a `SNES` evaluation function that would be passed to `SNESSetFunction()`
4529bcc50f1SBarry Smith 
4539bcc50f1SBarry Smith   Calling Sequence:
4549bcc50f1SBarry Smith + snes  - `SNES` context
4559bcc50f1SBarry Smith . u   - input vector
4569bcc50f1SBarry Smith . F   - function vector
4579bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4589bcc50f1SBarry Smith 
4599bcc50f1SBarry Smith   Level: beginner
4609bcc50f1SBarry Smith 
4618434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
4629bcc50f1SBarry Smith S*/
463*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESFunctionFn(SNES snes, Vec u, Vec F, void *ctx);
4649bcc50f1SBarry Smith 
4659bcc50f1SBarry Smith /*S
4668434afd1SBarry Smith   SNESObjectiveFn - A prototype of a `SNES` objective evaluation function that would be passed to `SNESSetObjective()`
4679bcc50f1SBarry Smith 
4689bcc50f1SBarry Smith   Calling Sequence:
4699bcc50f1SBarry Smith + snes  - `SNES` context
4709bcc50f1SBarry Smith . u   - input vector
4719bcc50f1SBarry Smith . o   - output value
4729bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4739bcc50f1SBarry Smith 
4749bcc50f1SBarry Smith   Level: beginner
4759bcc50f1SBarry Smith 
4768434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
4779bcc50f1SBarry Smith S*/
478*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESObjectiveFn(SNES snes, Vec u, PetscReal *o, void *ctx);
4799bcc50f1SBarry Smith 
4809bcc50f1SBarry Smith /*S
4818434afd1SBarry Smith   SNESJacobianFn - A prototype of a `SNES` Jacobian evaluation function that would be passed to `SNESSetJacobian()`
4829bcc50f1SBarry Smith 
4839bcc50f1SBarry Smith   Calling Sequence:
4849bcc50f1SBarry Smith + snes   - the `SNES` context obtained from `SNESCreate()`
4859bcc50f1SBarry Smith . u    - input vector
4869bcc50f1SBarry Smith . Amat - (approximate) Jacobian matrix
4879bcc50f1SBarry Smith . Pmat - matrix used to construct the preconditioner, often the same as `Amat`
4889bcc50f1SBarry Smith - ctx  - [optional] user-defined context for matrix evaluation routine
4899bcc50f1SBarry Smith 
4909bcc50f1SBarry Smith   Level: beginner
4919bcc50f1SBarry Smith 
4928434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESNGSFn`
4939bcc50f1SBarry Smith S*/
494*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESJacobianFn(SNES snes, Vec u, Mat Amat, Mat Pmat, void *ctx);
4959bcc50f1SBarry Smith 
4969bcc50f1SBarry Smith /*S
4978434afd1SBarry Smith   SNESNGSFn - A prototype of a `SNES` nonlinear Gauss-Seidel function that would be passed to `SNESSetNGS()`
4989bcc50f1SBarry Smith 
4999bcc50f1SBarry Smith   Calling Sequence:
5009bcc50f1SBarry Smith + snes   - the `SNES` context obtained from `SNESCreate()`
5019bcc50f1SBarry Smith . u    - the current solution, updated in place
502dd8e379bSPierre Jolivet . b    - the right-hand side vector (which may be `NULL`)
5039bcc50f1SBarry Smith - ctx  - [optional] user-defined context for matrix evaluation routine
5049bcc50f1SBarry Smith 
5059bcc50f1SBarry Smith   Level: beginner
5069bcc50f1SBarry Smith 
5078434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`
5089bcc50f1SBarry Smith S*/
509*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESNGSFn(SNES snes, Vec u, Vec b, void *ctx);
5109bcc50f1SBarry Smith 
51153e5d35bSStefano Zampini /*S
51253e5d35bSStefano Zampini   SNESUpdateFn - A prototype of a `SNES` update function that would be passed to `SNESSetUpdate()`
51353e5d35bSStefano Zampini 
51453e5d35bSStefano Zampini   Calling Sequence:
51553e5d35bSStefano Zampini + snes - `SNES` context
51653e5d35bSStefano Zampini - step - the current iteration index
51753e5d35bSStefano Zampini 
51853e5d35bSStefano Zampini   Level: advanced
51953e5d35bSStefano Zampini 
52053e5d35bSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESSetUpdate()`
52153e5d35bSStefano Zampini S*/
522*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESUpdateFn(SNES snes, PetscInt step);
52353e5d35bSStefano Zampini 
524b67197daSBarry Smith /* --------- Solving systems of nonlinear equations --------------- */
5258434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetFunction(SNES, Vec, SNESFunctionFn *, void *);
5268434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetFunction(SNES, Vec *, SNESFunctionFn **, void **);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESComputeFunction(SNES, Vec, Vec);
528bbc1464cSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeMFFunction(SNES, Vec, Vec);
52925acbd8eSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESSetInitialFunction(SNES, Vec);
53032f3f7c2SPeter Brune 
5318434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetJacobian(SNES, Mat, Mat, SNESJacobianFn *, void *);
5328434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetJacobian(SNES, Mat *, Mat *, SNESJacobianFn **, void **);
5338434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESObjectiveComputeFunctionDefaultFD;
5348434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefault;
5358434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefaultColor;
53678e7fe0eSHong Zhang PETSC_EXTERN PetscErrorCode SNESPruneJacobianColor(SNES, Mat, Mat);
5378434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetComputeInitialGuess(SNES, SNESInitialGuessFn *, void *);
5388434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetPicard(SNES, Vec, SNESFunctionFn *, Mat, Mat, SNESJacobianFn *, void *);
5398434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetPicard(SNES, Vec *, SNESFunctionFn **, Mat *, Mat *, SNESJacobianFn **, void **);
5408434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeFunction;
5418434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeMFFunction;
5428434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESPicardComputeJacobian;
54317bae607SBarry Smith 
5448434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetObjective(SNES, SNESObjectiveFn *, void *);
5458434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetObjective(SNES, SNESObjectiveFn **, void **);
5462a4ee8f2SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeObjective(SNES, Vec, PetscReal *);
5472a4ee8f2SPeter Brune 
54853e5d35bSStefano Zampini PETSC_EXTERN PetscErrorCode SNESSetUpdate(SNES, SNESUpdateFn *);
54953e5d35bSStefano Zampini 
550534ebe21SPeter Brune /*E
55116a05f60SBarry Smith    SNESNormSchedule - Frequency with which the norm is computed during a nonliner solve
552534ebe21SPeter Brune 
553dc4c0fb0SBarry Smith    Values:
554dc4c0fb0SBarry Smith +   `SNES_NORM_DEFAULT`            - use the default behavior for the current `SNESType`
555dc4c0fb0SBarry Smith .   `SNES_NORM_NONE`               - avoid all norm computations
556dc4c0fb0SBarry Smith .   `SNES_NORM_ALWAYS`             - compute the norms whenever possible
557dc4c0fb0SBarry Smith .   `SNES_NORM_INITIAL_ONLY`       - compute the norm only when the algorithm starts
558dc4c0fb0SBarry Smith .   `SNES_NORM_FINAL_ONLY`         - compute the norm only when the algorithm finishes
559dc4c0fb0SBarry Smith -   `SNES_NORM_INITIAL_FINAL_ONLY` - compute the norm at the start and end of the algorithm
560534ebe21SPeter Brune 
56116a05f60SBarry Smith    Level: advanced
56216a05f60SBarry Smith 
563534ebe21SPeter Brune    Notes:
564dc4c0fb0SBarry Smith    Support for these is highly dependent on the solver.
565dc4c0fb0SBarry Smith 
566dc4c0fb0SBarry Smith    Some options limit the convergence tests that can be used.
567dc4c0fb0SBarry Smith 
568dc4c0fb0SBarry Smith    The `SNES_NORM_NONE` option is most commonly used when the nonlinear solver is being used as a smoother, for example for `SNESFAS`
569dc4c0fb0SBarry Smith 
570534ebe21SPeter Brune    This is primarily used to turn off extra norm and function computation
571534ebe21SPeter Brune    when the solvers are composed.
572534ebe21SPeter Brune 
5731cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
574db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
575534ebe21SPeter Brune E*/
5769371c9d4SSatish Balay typedef enum {
5779371c9d4SSatish Balay   SNES_NORM_DEFAULT            = -1,
578fdacfa88SPeter Brune   SNES_NORM_NONE               = 0,
579365a6726SPeter Brune   SNES_NORM_ALWAYS             = 1,
580fdacfa88SPeter Brune   SNES_NORM_INITIAL_ONLY       = 2,
581fdacfa88SPeter Brune   SNES_NORM_FINAL_ONLY         = 3,
5829371c9d4SSatish Balay   SNES_NORM_INITIAL_FINAL_ONLY = 4
5839371c9d4SSatish Balay } SNESNormSchedule;
584365a6726SPeter Brune PETSC_EXTERN const char *const *const SNESNormSchedules;
5851957e957SBarry Smith 
586534ebe21SPeter Brune /*MC
58716a05f60SBarry Smith    SNES_NORM_NONE - Don't compute function and its L2 norm when possible
588534ebe21SPeter Brune 
589534ebe21SPeter Brune    Level: advanced
590534ebe21SPeter Brune 
591dc4c0fb0SBarry Smith    Note:
592534ebe21SPeter Brune    This is most useful for stationary solvers with a fixed number of iterations used as smoothers.
593534ebe21SPeter Brune 
5941cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_DEFAULT`
595534ebe21SPeter Brune M*/
596534ebe21SPeter Brune 
597fdacfa88SPeter Brune /*MC
598365a6726SPeter Brune    SNES_NORM_ALWAYS - Compute the function and its L2 norm at each iteration.
599534ebe21SPeter Brune 
600fdacfa88SPeter Brune    Level: advanced
601fdacfa88SPeter Brune 
602dc4c0fb0SBarry Smith    Note:
603fdacfa88SPeter Brune    Most solvers will use this no matter what norm type is passed to them.
604fdacfa88SPeter Brune 
6051cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_NONE`
606fdacfa88SPeter Brune M*/
607534ebe21SPeter Brune 
608534ebe21SPeter Brune /*MC
609534ebe21SPeter Brune    SNES_NORM_INITIAL_ONLY - Compute the function and its L2 at iteration 0, but do not update it.
610534ebe21SPeter Brune 
611534ebe21SPeter Brune    Level: advanced
612534ebe21SPeter Brune 
613534ebe21SPeter Brune    Notes:
614dc4c0fb0SBarry Smith    This method is useful in composed methods, when a true solution might actually be found before `SNESSolve()` is called.
615534ebe21SPeter Brune    This option enables the solve to abort on the zeroth iteration if this is the case.
616534ebe21SPeter Brune 
617534ebe21SPeter Brune    For solvers that require the computation of the L2 norm of the function as part of the method, this merely cancels
618534ebe21SPeter Brune    the norm computation at the last iteration (if possible).
619534ebe21SPeter Brune 
6201cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_FINAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
621534ebe21SPeter Brune M*/
622534ebe21SPeter Brune 
623534ebe21SPeter Brune /*MC
624534ebe21SPeter Brune    SNES_NORM_FINAL_ONLY - Compute the function and its L2 norm on only the final iteration.
625534ebe21SPeter Brune 
626534ebe21SPeter Brune    Level: advanced
627534ebe21SPeter Brune 
628dc4c0fb0SBarry Smith    Note:
629534ebe21SPeter Brune    For solvers that require the computation of the L2 norm of the function as part of the method, behaves
63087497f52SBarry Smith    exactly as `SNES_NORM_DEFAULT`.  This method is useful when the function is gotten after `SNESSolve()` and
631534ebe21SPeter Brune    used in subsequent computation for methods that do not need the norm computed during the rest of the
632534ebe21SPeter Brune    solution procedure.
633534ebe21SPeter Brune 
6341cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_INITIAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
635534ebe21SPeter Brune M*/
636534ebe21SPeter Brune 
637534ebe21SPeter Brune /*MC
638534ebe21SPeter Brune    SNES_NORM_INITIAL_FINAL_ONLY - Compute the function and its L2 norm on only the initial and final iterations.
639534ebe21SPeter Brune 
640534ebe21SPeter Brune    Level: advanced
641534ebe21SPeter Brune 
642dc4c0fb0SBarry Smith    Note:
64387497f52SBarry Smith    This method combines the benefits of `SNES_NORM_INITIAL_ONLY` and `SNES_NORM_FINAL_ONLY`.
644534ebe21SPeter Brune 
6451cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_SNES_NORM_INITIAL_ONLY`, `SNES_NORM_FINAL_ONLY`
646534ebe21SPeter Brune M*/
647534ebe21SPeter Brune 
648365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetNormSchedule(SNES, SNESNormSchedule);
649365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetNormSchedule(SNES, SNESNormSchedule *);
650c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetFunctionNorm(SNES, PetscReal);
651c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESGetFunctionNorm(SNES, PetscReal *);
652c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetUpdateNorm(SNES, PetscReal *);
653c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetSolutionNorm(SNES, PetscReal *);
654534ebe21SPeter Brune 
65547073ea2SPeter Brune /*E
65647073ea2SPeter Brune    SNESFunctionType - Type of function computed
65747073ea2SPeter Brune 
658dc4c0fb0SBarry Smith    Values:
659dc4c0fb0SBarry Smith +  `SNES_FUNCTION_DEFAULT`          - the default behavior for the current `SNESType`
660dc4c0fb0SBarry Smith .  `SNES_FUNCTION_UNPRECONDITIONED` - the original function provided
661dc4c0fb0SBarry Smith -  `SNES_FUNCTION_PRECONDITIONED`   - the modification of the function by the preconditioner
66247073ea2SPeter Brune 
66316a05f60SBarry Smith    Level: advanced
66416a05f60SBarry Smith 
665dc4c0fb0SBarry Smith    Note:
666dc4c0fb0SBarry Smith    Support for these is dependent on the solver.
667dc4c0fb0SBarry Smith 
6681cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
669db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
67047073ea2SPeter Brune E*/
6719371c9d4SSatish Balay typedef enum {
6729371c9d4SSatish Balay   SNES_FUNCTION_DEFAULT          = -1,
67347073ea2SPeter Brune   SNES_FUNCTION_UNPRECONDITIONED = 0,
6749371c9d4SSatish Balay   SNES_FUNCTION_PRECONDITIONED   = 1
6759371c9d4SSatish Balay } SNESFunctionType;
67647073ea2SPeter Brune PETSC_EXTERN const char *const *const SNESFunctionTypes;
67747073ea2SPeter Brune 
67847073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetFunctionType(SNES, SNESFunctionType);
67947073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetFunctionType(SNES, SNESFunctionType *);
680f1c6b773SPeter Brune 
6818434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNGS(SNES, SNESNGSFn *, void *);
6828434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNGS(SNES, SNESNGSFn **, void **);
683be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGS(SNES, Vec, Vec);
684b6266c6eSPeter Brune 
685be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetSweeps(SNES, PetscInt);
686be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetSweeps(SNES, PetscInt *);
687be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt);
688be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
689badc63e7SPeter Brune 
6904fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES, PetscBool);
6914fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES, PetscBool *);
6924fc747eaSLawrence Mitchell 
6933ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode SNESShellGetContext(SNES, void *);
694014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESShellSetContext(SNES, void *);
695014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESShellSetSolve(SNES, PetscErrorCode (*)(SNES, Vec));
696646217ecSPeter Brune 
697c5ae4b9aSBarry Smith /* --------- Routines specifically for line search methods --------------- */
698c5ae4b9aSBarry Smith 
699872b6db9SPeter Brune /*S
70087497f52SBarry Smith    SNESLineSearch - Abstract PETSc object that manages line-search operations for nonlinear solvers
7019e764e56SPeter Brune 
7029e764e56SPeter Brune    Level: beginner
7039e764e56SPeter Brune 
7041cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNES`
7059e764e56SPeter Brune S*/
706907376e6SBarry Smith typedef struct _p_LineSearch *SNESLineSearch;
7079e764e56SPeter Brune 
7089e764e56SPeter Brune /*J
7090b4b7b1cSBarry Smith    SNESLineSearchType - String with the name of a PETSc line search method `SNESLineSearch`. Provides all the linesearches for the nonlinear solvers, `SNES`,
7100b4b7b1cSBarry Smith                         in PETSc.
7119e764e56SPeter Brune 
712ceaaa498SBarry Smith    Values:
713ceaaa498SBarry Smith +  `SNESLINESEARCHBASIC`   - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
714ceaaa498SBarry Smith .  `SNESLINESEARCHBT`      - Backtracking line search over the L2 norm of the function
715ceaaa498SBarry Smith .  `SNESLINESEARCHL2`      - Secant line search over the L2 norm of the function
7160b4b7b1cSBarry Smith .  `SNESLINESEARCHCP`      - Critical point secant line search assuming $F(x) = \nabla G(x)$ for some unknown $G(x)$
717ceaaa498SBarry Smith .  `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch
718ceaaa498SBarry Smith -  `SNESLINESEARCHSHELL`   - User provided `SNESLineSearch` implementation
719ceaaa498SBarry Smith 
72095bd0b28SBarry Smith    Level: beginner
72195bd0b28SBarry Smith 
7220b4b7b1cSBarry Smith    Note:
7230b4b7b1cSBarry Smith    Use `SNESLineSearchSetType()` or the options database key `-snes_linesearch_type` to set
7240b4b7b1cSBarry Smith    the specific line search algorithm to use with a given `SNES` object. Not all `SNESType` can utilize a line search.
7250b4b7b1cSBarry Smith 
7261cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetType()`, `SNES`
7279e764e56SPeter Brune J*/
72819fd82e9SBarry Smith typedef const char *SNESLineSearchType;
7299e764e56SPeter Brune #define SNESLINESEARCHBT        "bt"
730d4c6564cSPatrick Farrell #define SNESLINESEARCHNLEQERR   "nleqerr"
731c87759e9SPeter Brune #define SNESLINESEARCHBASIC     "basic"
7320b00b554SBarry Smith #define SNESLINESEARCHNONE      "none"
7339e764e56SPeter Brune #define SNESLINESEARCHL2        "l2"
7349e764e56SPeter Brune #define SNESLINESEARCHCP        "cp"
7359e764e56SPeter Brune #define SNESLINESEARCHSHELL     "shell"
736b5badacbSBarry Smith #define SNESLINESEARCHNCGLINEAR "ncglinear"
737b9d635d7SJonas Heinzmann #define SNESLINESEARCHBISECTION "bisection"
7389e764e56SPeter Brune 
739140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESList;
740014dd563SJed Brown PETSC_EXTERN PetscClassId      SNESLINESEARCH_CLASSID;
741140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESLineSearchList;
7429e764e56SPeter Brune 
743b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_LINEAR    1
744b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_QUADRATIC 2
745b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_CUBIC     3
7469e764e56SPeter Brune 
7479bcc50f1SBarry Smith /*S
7488434afd1SBarry Smith   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that projects a vector onto the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
7499bcc50f1SBarry Smith 
7509bcc50f1SBarry Smith   Calling Sequence:
7519bcc50f1SBarry Smith + snes  - `SNES` context
7529bcc50f1SBarry Smith - u     - the vector to project to the bounds
7539bcc50f1SBarry Smith 
7549bcc50f1SBarry Smith   Level: advanced
7559bcc50f1SBarry Smith 
7569bcc50f1SBarry Smith   Note:
7578434afd1SBarry Smith   The deprecated `SNESLineSearchVIProjectFunc` still works as a replacement for `SNESLineSearchVIProjectFn` *.
7589bcc50f1SBarry Smith 
7599bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
7609bcc50f1SBarry Smith S*/
761*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode             SNESLineSearchVIProjectFn(SNES snes, Vec u);
762*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVIProjectFn *SNESLineSearchVIProjectFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVIProjectFn*", );
7639bcc50f1SBarry Smith 
7649bcc50f1SBarry Smith /*S
7658434afd1SBarry Smith   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that computes the norm of the active set variables in a vector in a VI solve,
7669bcc50f1SBarry Smith   passed to `SNESLineSearchSetVIFunctions()`
7679bcc50f1SBarry Smith 
7689bcc50f1SBarry Smith   Calling Sequence:
7699bcc50f1SBarry Smith + snes  - `SNES` context
7709bcc50f1SBarry Smith . f     - the vector to compute the norm of
7719bcc50f1SBarry Smith . u     - the current solution, entries that are on the VI bounds are ignored
7729bcc50f1SBarry Smith - fnorm - the resulting norm
7739bcc50f1SBarry Smith 
7749bcc50f1SBarry Smith   Level: advanced
7759bcc50f1SBarry Smith 
7769bcc50f1SBarry Smith   Note:
7778434afd1SBarry Smith   The deprecated `SNESLineSearchVINormFunc` still works as a replacement for `SNESLineSearchVINormFn` *.
7789bcc50f1SBarry Smith 
7799bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
7809bcc50f1SBarry Smith S*/
781*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode          SNESLineSearchVINormFn(SNES snes, Vec f, Vec u, PetscReal *fnorm);
782*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVINormFn *SNESLineSearchVINormFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVINormFnn*", );
7839bcc50f1SBarry Smith 
784d5def619SJonas Heinzmann /*S
785d5def619SJonas Heinzmann   SNESLineSearchVIDirDerivFn - A prototype of a `SNES` function that computes the directional derivative considering the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
786d5def619SJonas Heinzmann 
787d5def619SJonas Heinzmann   Calling Sequence:
788d5def619SJonas Heinzmann + snes  - `SNES` context
789d5def619SJonas Heinzmann . f     - the function vector to compute the directional derivative with
790d5def619SJonas Heinzmann . u     - the current solution, entries that are on the VI bounds are ignored
791d5def619SJonas Heinzmann . y     - the direction to compute the directional derivative
792d5def619SJonas Heinzmann - fty   - the resulting directional derivative
793d5def619SJonas Heinzmann 
794d5def619SJonas Heinzmann   Level: advanced
795d5def619SJonas Heinzmann 
796d5def619SJonas Heinzmann .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchVIProjectFn`, `SNESLineSearchVIProjectFn`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetVIFunctions()`
797d5def619SJonas Heinzmann S*/
798*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchVIDirDerivFn(SNES snes, Vec f, Vec u, Vec y, PetscScalar *fty);
799d5def619SJonas Heinzmann 
800*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchApplyFn(SNESLineSearch);
801*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchApplyFn      *SNESLineSearchApplyFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
802*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode              SNESLineSearchShellApplyFn(SNESLineSearch, void *);
803*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchShellApplyFn *SNESLineSearchUserFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
8049e764e56SPeter Brune 
805014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate(MPI_Comm, SNESLineSearch *);
806014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchReset(SNESLineSearch);
807014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchView(SNESLineSearch, PetscViewer);
808014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *);
809a80ff896SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetType(SNESLineSearch, SNESLineSearchType *);
81019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetType(SNESLineSearch, SNESLineSearchType);
811014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch);
812ed07d7d7SPeter Brune PETSC_EXTERN PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch, PetscErrorCode (*)(SNES, Vec, Vec));
813014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetUp(SNESLineSearch);
814014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchApply(SNESLineSearch, Vec, Vec, PetscReal *, Vec);
815014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch, Vec, Vec, PetscBool *);
816014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *);
817fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch, PetscInt);
8189e764e56SPeter Brune 
8199bd66eb0SPeter Brune /* set the functions for precheck and postcheck */
82086d74e61SPeter Brune 
821f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx);
822f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);
8239e764e56SPeter Brune 
824f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx);
825f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);
8269e764e56SPeter Brune 
8279bd66eb0SPeter Brune /* set the functions for VI-specific line search operations */
8289bd66eb0SPeter Brune 
829d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn *, SNESLineSearchVINormFn *, SNESLineSearchVIDirDerivFn *);
830d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn **, SNESLineSearchVINormFn **, SNESLineSearchVIDirDerivFn **);
8319bd66eb0SPeter Brune 
8329e764e56SPeter Brune /* pointers to the associated SNES in order to be able to get the function evaluation out */
833014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch, SNES);
834014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch, SNES *);
8359e764e56SPeter Brune 
8369e764e56SPeter Brune /* set and get the parameters and vectors */
837014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
838014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscInt);
8399e764e56SPeter Brune 
840014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch, Vec, Vec, PetscBool *, void *);
84186d74e61SPeter Brune 
842014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch, PetscReal *);
843014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch, PetscReal);
8449e764e56SPeter Brune 
845014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch, PetscReal *);
846014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch, PetscReal);
8479e764e56SPeter Brune 
848ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch, PetscInt *);
849ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch, PetscInt);
85059405d5eSPeter Brune 
851422a814eSBarry Smith /*E
852dc4c0fb0SBarry Smith     SNESLineSearchReason - indication if the line search has succeeded or failed and why
853422a814eSBarry Smith 
854dc4c0fb0SBarry Smith   Values:
855dc4c0fb0SBarry Smith +  `SNES_LINESEARCH_SUCCEEDED`       - the line search succeeded
856dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_NANORINF` - a not a number of infinity appeared in the computions
857dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_DOMAIN`   - the function was evaluated outside of its domain, see `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
858dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_REDUCT`   - the linear search failed to get the requested decrease in its norm or objective
859dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_USER`     - used by `SNESLINESEARCHNLEQERR` to indicate the user changed the search direction inappropriately
860dc4c0fb0SBarry Smith -  `SNES_LINESEARCH_FAILED_FUNCTION` - indicates the maximum number of function evaluations allowed has been surpassed, `SNESConvergedReason` is also
861dc4c0fb0SBarry Smith                                        set to `SNES_DIVERGED_FUNCTION_COUNT`
862422a814eSBarry Smith 
86316a05f60SBarry Smith    Level: intermediate
86416a05f60SBarry Smith 
865dc4c0fb0SBarry Smith    Developer Note:
866dc4c0fb0SBarry Smith    Some of these reasons overlap with values of `SNESConvergedReason`
867422a814eSBarry Smith 
8681cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`,
869dc4c0fb0SBarry Smith           `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
870422a814eSBarry Smith E*/
8719371c9d4SSatish Balay typedef enum {
8729371c9d4SSatish Balay   SNES_LINESEARCH_SUCCEEDED,
873422a814eSBarry Smith   SNES_LINESEARCH_FAILED_NANORINF,
874e9b602ebSSatish Balay   SNES_LINESEARCH_FAILED_DOMAIN,
875d5b43468SJose E. Roman   SNES_LINESEARCH_FAILED_REDUCT, /* INSUFFICIENT REDUCTION */
876e9b602ebSSatish Balay   SNES_LINESEARCH_FAILED_USER,
8779371c9d4SSatish Balay   SNES_LINESEARCH_FAILED_FUNCTION
8789371c9d4SSatish Balay } SNESLineSearchReason;
879422a814eSBarry Smith 
880422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetReason(SNESLineSearch, SNESLineSearchReason *);
881422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetReason(SNESLineSearch, SNESLineSearchReason);
8829e764e56SPeter Brune 
883014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch, Vec *, Vec *, Vec *, Vec *, Vec *);
884014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch, Vec, Vec, Vec, Vec, Vec);
8859e764e56SPeter Brune 
886014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *);
887014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch, PetscReal, PetscReal, PetscReal);
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch, PetscBool);
8909e764e56SPeter Brune 
891dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitor(SNESLineSearch);
89249abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, void *), void *, PetscCtxDestroyFn *);
893d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch, const char[], const char[], const char[], PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *));
894dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch);
895dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch, PetscViewer);
896dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch, PetscViewer *);
897d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch, PetscViewerAndFormat *);
8989e764e56SPeter Brune 
899ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch, const char[]);
900ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch, const char *[]);
9019e764e56SPeter Brune 
9029e764e56SPeter Brune /* Shell interface functions */
903*21bb954dSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellSetApply(SNESLineSearch, SNESLineSearchShellApplyFn *, void *);
9048434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellGetApply(SNESLineSearch, SNESLineSearchShellApplyFn **, void **);
9059bcc50f1SBarry Smith 
906*21bb954dSBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellSetApply()", ) static inline PetscErrorCode SNESLineSearchShellSetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn *f, void *ctx)
9079bcc50f1SBarry Smith {
9089bcc50f1SBarry Smith   return SNESLineSearchShellSetApply(ls, f, ctx);
9099bcc50f1SBarry Smith }
9109bcc50f1SBarry Smith 
911*21bb954dSBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellGetApply()", ) static inline PetscErrorCode SNESLineSearchShellGetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn **f, void **ctx)
9129bcc50f1SBarry Smith {
9139bcc50f1SBarry Smith   return SNESLineSearchShellGetApply(ls, f, ctx);
9149bcc50f1SBarry Smith }
9159e764e56SPeter Brune 
9162f4102e2SPeter Brune /* BT interface functions */
917014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch, PetscReal);
918014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch, PetscReal *);
9192f4102e2SPeter Brune 
9209e764e56SPeter Brune /*register line search types */
921bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchRegister(const char[], PetscErrorCode (*)(SNESLineSearch));
9229e764e56SPeter Brune 
923720c9a41SShri Abhyankar /* Routines for VI solver */
924014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetVariableBounds(SNES, Vec, Vec);
925cf836535SBlaise Bourdin PETSC_EXTERN PetscErrorCode SNESVIGetVariableBounds(SNES, Vec *, Vec *);
926014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
927014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetInactiveSet(SNES, IS *);
928014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetActiveSetIS(SNES, Vec, Vec, IS *);
929014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES, Vec, Vec, PetscReal *);
930d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFtY(SNES, Vec, Vec, Vec, PetscScalar *);
931014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetRedundancyCheck(SNES, PetscErrorCode (*)(SNES, IS, IS *, void *), void *);
932f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeMeritFunction(Vec, PetscReal *, PetscReal *);
933f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeFunction(SNES, Vec, Vec, void *);
934f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMSetVI(DM, IS);
935f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMDestroyVI(DM);
936f5ea5bd2SShri Abhyankar 
937014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTestLocalMin(SNES);
938da9b6338SBarry Smith 
939eef9c623SLois Curfman McInnes /* Should this routine be private? */
940d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian(SNES, Vec, Mat, Mat);
941e885f1abSBarry Smith PETSC_EXTERN PetscErrorCode SNESTestJacobian(SNES);
942494a190aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESTestFunction(SNES);
943ddbbdb52SLois Curfman McInnes 
944014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetDM(SNES, DM);
945014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetDM(SNES, DM *);
946be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPC(SNES, SNES);
947be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPC(SNES, SNES *);
9483ad1a0b9SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESHasNPC(SNES, PetscBool *);
949be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESApplyNPC(SNES, Vec, Vec, Vec);
950be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCFunction(SNES, Vec, PetscReal *);
951be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeFunctionDefaultNPC(SNES, Vec, Vec);
952be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPCSide(SNES, PCSide);
953be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCSide(SNES, PCSide *);
9547601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLineSearch(SNES, SNESLineSearch);
9557601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLineSearch(SNES, SNESLineSearch *);
9566c699258SBarry Smith 
957edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESGetLineSearch()", ) static inline PetscErrorCode SNESGetSNESLineSearch(SNES snes, SNESLineSearch *ls)
958d71ae5a4SJacob Faibussowitsch {
9599371c9d4SSatish Balay   return SNESGetLineSearch(snes, ls);
9609371c9d4SSatish Balay }
961edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESSetLineSearch()", ) static inline PetscErrorCode SNESSetSNESLineSearch(SNES snes, SNESLineSearch ls)
962d71ae5a4SJacob Faibussowitsch {
9639371c9d4SSatish Balay   return SNESSetLineSearch(snes, ls);
9649371c9d4SSatish Balay }
9658b7b3213SJed Brown 
966014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetUpMatrices(SNES);
9678434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunction(DM, SNESFunctionFn *, void *);
9688434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetFunction(DM, SNESFunctionFn **, void **);
96949abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionContextDestroy(DM, PetscCtxDestroyFn *);
9708434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetMFFunction(DM, SNESFunctionFn *, void *);
9718434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetNGS(DM, SNESNGSFn *, void *);
9728434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetNGS(DM, SNESNGSFn **, void **);
9738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobian(DM, SNESJacobianFn *, void *);
9748434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetJacobian(DM, SNESJacobianFn **, void **);
97549abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianContextDestroy(DM, PetscCtxDestroyFn *);
9768434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetPicard(DM, SNESFunctionFn *, SNESJacobianFn *, void *);
9778434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetPicard(DM, SNESFunctionFn **, SNESJacobianFn **, void **);
9788434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetObjective(DM, SNESObjectiveFn *, void *);
9798434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetObjective(DM, SNESObjectiveFn **, void **);
9804dd50a75SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMSNES(DM, DM);
9816cab3a1bSJed Brown 
982*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionFn(DMDALocalInfo *, void *, void *, void *);
983*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianFn(DMDALocalInfo *, void *, Mat, Mat, void *);
984*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveFn(DMDALocalInfo *, void *, PetscReal *, void *);
985c9d099b5SPeter Brune 
986*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionVecFn(DMDALocalInfo *, Vec, Vec, void *);
987*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianVecFn(DMDALocalInfo *, Vec, Mat, Mat, void *);
988*21bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveVecFn(DMDALocalInfo *, Vec, PetscReal *, void *);
9895505f7afSJunchao Zhang 
9908434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocal(DM, InsertMode, DMDASNESFunctionFn *, void *);
9918434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocal(DM, DMDASNESJacobianFn *, void *);
9928434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocal(DM, DMDASNESObjectiveFn *, void *);
9938434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetPicardLocal(DM, InsertMode, DMDASNESFunctionFn *, DMDASNESJacobianFn, void *);
9946cab3a1bSJed Brown 
9958434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocalVec(DM, InsertMode, DMDASNESFunctionVecFn *, void *);
9968434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocalVec(DM, DMDASNESJacobianVecFn *, void *);
9978434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocalVec(DM, DMDASNESObjectiveVecFn *, void *);
9985505f7afSJunchao Zhang 
999bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMSNESSetBoundaryLocal(DM, PetscErrorCode (*)(DM, Vec, void *), void *);
10006493148fSStefano Zampini PETSC_EXTERN PetscErrorCode DMSNESSetObjectiveLocal(DM, PetscErrorCode (*)(DM, Vec, PetscReal *, void *), void *);
1001ff35dfedSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionLocal(DM, PetscErrorCode (*)(DM, Vec, Vec, void *), void *);
1002d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianLocal(DM, PetscErrorCode (*)(DM, Vec, Mat, Mat, void *), void *);
100328d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetBoundaryLocal(DM, PetscErrorCode (**)(DM, Vec, void *), void **);
10046493148fSStefano Zampini PETSC_EXTERN PetscErrorCode DMSNESGetObjectiveLocal(DM, PetscErrorCode (**)(DM, Vec, PetscReal *, void *), void **);
100528d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetFunctionLocal(DM, PetscErrorCode (**)(DM, Vec, Vec, void *), void **);
100628d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetJacobianLocal(DM, PetscErrorCode (**)(DM, Vec, Mat, Mat, void *), void **);
1007ff35dfedSBarry Smith 
1008b665c654SMatthew G Knepley /* Routines for Multiblock solver */
1009014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetFields(SNES, const char[], PetscInt, const PetscInt *);
1010014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetIS(SNES, const char[], IS);
1011014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
1012014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);
10134bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMultiblockGetSubSNES(SNES, PetscInt *, SNES *[]);
1014aeed3662SMatthew G Knepley 
101537e1895aSJed Brown /*J
101687497f52SBarry Smith    SNESMSType - String with the name of a PETSc `SNESMS` method.
101737e1895aSJed Brown 
101837e1895aSJed Brown    Level: intermediate
101937e1895aSJed Brown 
10201cc06b55SBarry Smith .seealso: [](ch_snes), `SNESMS`, `SNESMSGetType()`, `SNESMSSetType()`, `SNES`
102137e1895aSJed Brown J*/
102219fd82e9SBarry Smith typedef const char *SNESMSType;
102337e1895aSJed Brown #define SNESMSM62       "m62"
102437e1895aSJed Brown #define SNESMSEULER     "euler"
1025a97cb6bcSJed Brown #define SNESMSJAMESON83 "jameson83"
10263847c725SLisandro Dalcin #define SNESMSVLTP11    "vltp11"
1027a97cb6bcSJed Brown #define SNESMSVLTP21    "vltp21"
1028a97cb6bcSJed Brown #define SNESMSVLTP31    "vltp31"
1029a97cb6bcSJed Brown #define SNESMSVLTP41    "vltp41"
1030a97cb6bcSJed Brown #define SNESMSVLTP51    "vltp51"
1031a97cb6bcSJed Brown #define SNESMSVLTP61    "vltp61"
103237e1895aSJed Brown 
103319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSRegister(SNESMSType, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[]);
10344bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMSRegisterAll(void);
103557715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetType(SNES, SNESMSType *);
103619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSSetType(SNES, SNESMSType);
103757715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetDamping(SNES, PetscReal *);
103857715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSSetDamping(SNES, PetscReal);
1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSFinalizePackage(void);
1040607a6623SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSInitializePackage(void);
1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSRegisterDestroy(void);
104237e1895aSJed Brown 
1043ceaaa498SBarry Smith /*MC
1044ceaaa498SBarry Smith    SNESNGMRESRestartType - the restart approach used by `SNESNGMRES`
104513a62661SPeter Brune 
1046ceaaa498SBarry Smith   Values:
1047ceaaa498SBarry Smith +   `SNES_NGMRES_RESTART_NONE`       - never restart
1048ceaaa498SBarry Smith .   `SNES_NGMRES_RESTART_DIFFERENCE` - restart based upon difference criteria
1049ceaaa498SBarry Smith -   `SNES_NGMRES_RESTART_PERIODIC`   - restart after a fixed number of iterations
1050ceaaa498SBarry Smith 
1051ceaaa498SBarry Smith   Options Database Keys:
1052ceaaa498SBarry Smith + -snes_ngmres_restart_type <difference,periodic,none> - set the restart type
1053af27ebaaSBarry Smith - -snes_ngmres_restart <30>                            - sets the number of iterations before restart for periodic
1054ceaaa498SBarry Smith 
105595bd0b28SBarry Smith    Level: intermediate
105695bd0b28SBarry Smith 
1057ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1058ceaaa498SBarry Smith           `SNESNGMRESGetRestartType()`, `SNESNGMRESSelectType`
1059ceaaa498SBarry Smith M*/
106013a62661SPeter Brune typedef enum {
106113a62661SPeter Brune   SNES_NGMRES_RESTART_NONE       = 0,
106213a62661SPeter Brune   SNES_NGMRES_RESTART_PERIODIC   = 1,
10639371c9d4SSatish Balay   SNES_NGMRES_RESTART_DIFFERENCE = 2
10649371c9d4SSatish Balay } SNESNGMRESRestartType;
10656a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];
106613a62661SPeter Brune 
1067ceaaa498SBarry Smith /*MC
1068ceaaa498SBarry Smith    SNESNGMRESSelectType - the approach used by `SNESNGMRES` to determine how the candidate solution and
1069ceaaa498SBarry Smith   combined solution are used to create the next iterate.
1070ceaaa498SBarry Smith 
1071ceaaa498SBarry Smith    Values:
1072ceaaa498SBarry Smith +   `SNES_NGMRES_SELECT_NONE`       - choose the combined solution all the time
1073ceaaa498SBarry Smith .   `SNES_NGMRES_SELECT_DIFFERENCE` - choose based upon the selection criteria
1074ceaaa498SBarry Smith -   `SNES_NGMRES_SELECT_LINESEARCH` - choose based upon line search combination
1075ceaaa498SBarry Smith 
1076ceaaa498SBarry Smith   Options Database Key:
1077ceaaa498SBarry Smith . -snes_ngmres_select_type<difference,none,linesearch> - select type
1078ceaaa498SBarry Smith 
107995bd0b28SBarry Smith    Level: intermediate
108095bd0b28SBarry Smith 
1081ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1082ceaaa498SBarry Smith           `SNESNGMRESGetRestartType()`, `SNESNGMRESRestartType`
1083ceaaa498SBarry Smith M*/
108413a62661SPeter Brune typedef enum {
108513a62661SPeter Brune   SNES_NGMRES_SELECT_NONE       = 0,
108613a62661SPeter Brune   SNES_NGMRES_SELECT_DIFFERENCE = 1,
10879371c9d4SSatish Balay   SNES_NGMRES_SELECT_LINESEARCH = 2
10889371c9d4SSatish Balay } SNESNGMRESSelectType;
10896a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESSelectTypes[];
109013a62661SPeter Brune 
1091014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartType(SNES, SNESNGMRESRestartType);
1092014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetSelectType(SNES, SNESNGMRESSelectType);
109323b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartFmRise(SNES, PetscBool);
109423b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESGetRestartFmRise(SNES, PetscBool *);
109513a62661SPeter Brune 
1096ceaaa498SBarry Smith /*MC
1097ceaaa498SBarry Smith    SNESNCGType - the conjugate update approach for `SNESNCG`
10980a844d1aSPeter Brune 
1099ceaaa498SBarry Smith    Values:
1100ceaaa498SBarry Smith +   `SNES_NCG_FR`  - Fletcher-Reeves update
1101ceaaa498SBarry Smith .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update, the default and the only one that tolerates generalized search directions
1102ceaaa498SBarry Smith .   `SNES_NCG_HS`  - Hestenes-Steifel update
1103ceaaa498SBarry Smith .   `SNES_NCG_DY`  - Dai-Yuan update
1104ceaaa498SBarry Smith -   `SNES_NCG_CD`  - Conjugate Descent update
1105ceaaa498SBarry Smith 
1106ceaaa498SBarry Smith   Options Database Key:
1107ceaaa498SBarry Smith . -snes_ncg_type<fr,prp,hs,dy,cd> - select type
1108ceaaa498SBarry Smith 
110995bd0b28SBarry Smith    Level: intermediate
111095bd0b28SBarry Smith 
1111ceaaa498SBarry Smith .seealso: `SNES, `SNESNCG`, `SNESNCGSetType()`
1112ceaaa498SBarry Smith M*/
11130a844d1aSPeter Brune typedef enum {
11140de8b71eSPeter Brune   SNES_NCG_FR  = 0,
11150de8b71eSPeter Brune   SNES_NCG_PRP = 1,
11160de8b71eSPeter Brune   SNES_NCG_HS  = 2,
11170de8b71eSPeter Brune   SNES_NCG_DY  = 3,
11189371c9d4SSatish Balay   SNES_NCG_CD  = 4
11199371c9d4SSatish Balay } SNESNCGType;
11206a6fc655SJed Brown PETSC_EXTERN const char *const SNESNCGTypes[];
11210a844d1aSPeter Brune 
1122014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNCGSetType(SNES, SNESNCGType);
112313a62661SPeter Brune 
1124ceaaa498SBarry Smith /*MC
1125ceaaa498SBarry Smith    SNESQNScaleType - the scaling type used by `SNESQN`
1126ceaaa498SBarry Smith 
1127ceaaa498SBarry Smith    Values:
1128ceaaa498SBarry Smith +   `SNES_QN_SCALE_NONE`     - don't scale the problem
1129ceaaa498SBarry Smith .   `SNES_QN_SCALE_SCALAR`   - use Shanno scaling
1130ceaaa498SBarry Smith .   `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
1131ceaaa498SBarry Smith -   `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with `SNESSetJacobian()`
1132ceaaa498SBarry Smith                                computed at the first iteration of `SNESQN` and at ever restart.
1133ceaaa498SBarry Smith 
1134ceaaa498SBarry Smith     Options Database Key:
1135ceaaa498SBarry Smith . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
1136ceaaa498SBarry Smith 
113795bd0b28SBarry Smith    Level: intermediate
113895bd0b28SBarry Smith 
1139ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNRestartType`
1140ceaaa498SBarry Smith M*/
11419371c9d4SSatish Balay typedef enum {
11429371c9d4SSatish Balay   SNES_QN_SCALE_DEFAULT  = 0,
11431efc8c45SPeter Brune   SNES_QN_SCALE_NONE     = 1,
114492f76d53SAlp Dener   SNES_QN_SCALE_SCALAR   = 2,
114592f76d53SAlp Dener   SNES_QN_SCALE_DIAGONAL = 3,
11469371c9d4SSatish Balay   SNES_QN_SCALE_JACOBIAN = 4
11479371c9d4SSatish Balay } SNESQNScaleType;
11486a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNScaleTypes[];
1149ceaaa498SBarry Smith 
1150ceaaa498SBarry Smith /*MC
1151ceaaa498SBarry Smith    SNESQNRestartType - the restart approached used by `SNESQN`
1152ceaaa498SBarry Smith 
1153ceaaa498SBarry Smith    Values:
1154ceaaa498SBarry Smith +   `SNES_QN_RESTART_NONE`     - never restart
1155ceaaa498SBarry Smith .   `SNES_QN_RESTART_POWELL`   - restart based upon descent criteria
1156ceaaa498SBarry Smith -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations
1157ceaaa498SBarry Smith 
1158ceaaa498SBarry Smith   Options Database Keys:
1159ceaaa498SBarry Smith + -snes_qn_restart_type <powell,periodic,none> - set the restart type
1160ceaaa498SBarry Smith - -snes_qn_m <m>                               - sets the number of stored updates and the restart period for periodic
1161ceaaa498SBarry Smith 
116295bd0b28SBarry Smith    Level: intermediate
116395bd0b28SBarry Smith 
1164ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNScaleType`
1165ceaaa498SBarry Smith M*/
11669371c9d4SSatish Balay typedef enum {
11679371c9d4SSatish Balay   SNES_QN_RESTART_DEFAULT  = 0,
11681efc8c45SPeter Brune   SNES_QN_RESTART_NONE     = 1,
11691efc8c45SPeter Brune   SNES_QN_RESTART_POWELL   = 2,
11709371c9d4SSatish Balay   SNES_QN_RESTART_PERIODIC = 3
11719371c9d4SSatish Balay } SNESQNRestartType;
11726a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNRestartTypes[];
1173ceaaa498SBarry Smith 
1174ceaaa498SBarry Smith /*MC
1175ceaaa498SBarry Smith    SNESQNType - the type used by `SNESQN`
1176ceaaa498SBarry Smith 
1177ceaaa498SBarry Smith   Values:
1178ceaaa498SBarry Smith +   `SNES_QN_LBFGS`      - LBFGS variant
1179ceaaa498SBarry Smith .   `SNES_QN_BROYDEN`    - Broyden variant
1180ceaaa498SBarry Smith -   `SNES_QN_BADBROYDEN` - Bad Broyden variant
1181ceaaa498SBarry Smith 
1182ceaaa498SBarry Smith   Options Database Key:
1183ceaaa498SBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
1184ceaaa498SBarry Smith 
118595bd0b28SBarry Smith    Level: intermediate
118695bd0b28SBarry Smith 
1187ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNSetType()`, `SNESQNScaleType`, `SNESQNRestartType`, `SNESQNSetRestartType()`
1188ceaaa498SBarry Smith M*/
11899371c9d4SSatish Balay typedef enum {
11909371c9d4SSatish Balay   SNES_QN_LBFGS      = 0,
11911efc8c45SPeter Brune   SNES_QN_BROYDEN    = 1,
11921efc8c45SPeter Brune   SNES_QN_BADBROYDEN = 2
11931efc8c45SPeter Brune } SNESQNType;
11941efc8c45SPeter Brune PETSC_EXTERN const char *const SNESQNTypes[];
11950c777b0cSPeter Brune 
11961efc8c45SPeter Brune PETSC_EXTERN PetscErrorCode SNESQNSetType(SNES, SNESQNType);
1197014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetScaleType(SNES, SNESQNScaleType);
1198014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetRestartType(SNES, SNESQNRestartType);
11990c777b0cSPeter Brune 
1200e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetType(SNES, PCASMType *);
1201e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetType(SNES, PCASMType);
1202ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomains(SNES, PetscInt *, SNES *[], VecScatter *[], VecScatter *[], VecScatter *[]);
1203ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMSetSubdomains(SNES, PetscInt, SNES[], VecScatter[], VecScatter[], VecScatter[]);
1204610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetDamping(SNES, PetscReal);
1205610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetDamping(SNES, PetscReal *);
1206ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomainVecs(SNES, PetscInt *, Vec *[], Vec *[], Vec *[], Vec *[]);
1207d728fb7dSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES, PetscBool);
1208d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetSNES(SNES, PetscInt, SNES *);
1209d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetNumber(SNES, PetscInt *);
1210f10b3e88SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMSetWeight(SNES, Vec);
12110c777b0cSPeter Brune 
1212ad05f2b7SBarry Smith /*E
1213ad05f2b7SBarry Smith   SNESCompositeType - Determines how two or more preconditioners are composed with the `SNESType` of `SNESCOMPOSITE`
1214ad05f2b7SBarry Smith 
1215ad05f2b7SBarry Smith   Values:
1216ad05f2b7SBarry Smith + `SNES_COMPOSITE_ADDITIVE`        - results from application of all preconditioners are added together
1217ad05f2b7SBarry Smith . `SNES_COMPOSITE_MULTIPLICATIVE`  - preconditioners are applied sequentially to the residual freshly
1218ad05f2b7SBarry Smith                                      computed after the previous preconditioner application
1219ad05f2b7SBarry Smith - `SNES_COMPOSITE_ADDITIVEOPTIMAL` - uses a linear combination of the solutions obtained with each preconditioner that approximately minimize the function
1220ad05f2b7SBarry Smith                                      value at the new iteration.
1221ad05f2b7SBarry Smith 
1222ad05f2b7SBarry Smith    Level: beginner
1223ad05f2b7SBarry Smith 
1224ad05f2b7SBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `PCCompositeType`
1225ad05f2b7SBarry Smith E*/
12269371c9d4SSatish Balay typedef enum {
12279371c9d4SSatish Balay   SNES_COMPOSITE_ADDITIVE,
12289371c9d4SSatish Balay   SNES_COMPOSITE_MULTIPLICATIVE,
12299371c9d4SSatish Balay   SNES_COMPOSITE_ADDITIVEOPTIMAL
12309371c9d4SSatish Balay } SNESCompositeType;
1231eed5f15bSPeter Brune PETSC_EXTERN const char *const SNESCompositeTypes[];
1232eed5f15bSPeter Brune 
1233eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetType(SNES, SNESCompositeType);
1234eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeAddSNES(SNES, SNESType);
1235eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeGetSNES(SNES, PetscInt, SNES *);
1236a6b47ab3SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESCompositeGetNumber(SNES, PetscInt *);
12378f626970SPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetDamping(SNES, PetscInt, PetscReal);
1238eed5f15bSPeter Brune 
123968e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetDiscretisationInfo(SNES, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
12404d04e9f1SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetComputeOperator(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
124139fd2e8aSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetComputeFunction(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
124268e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetConstructType(SNES, PCPatchConstructType, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *);
124368e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetCellNumbering(SNES, PetscSection);
1244561742edSMatthew G. Knepley 
1245a055c8caSBarry Smith /*E
1246a055c8caSBarry Smith     SNESFASType - Determines the type of nonlinear multigrid method that is run.
1247a055c8caSBarry Smith 
1248a055c8caSBarry Smith    Values:
124987497f52SBarry Smith +  `SNES_FAS_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `SNESFASSetCycles()`
125087497f52SBarry Smith .  `SNES_FAS_ADDITIVE`                 - additive FAS cycle
125187497f52SBarry Smith .  `SNES_FAS_FULL`                     - full FAS cycle
125287497f52SBarry Smith -  `SNES_FAS_KASKADE`                  - Kaskade FAS cycle
1253a055c8caSBarry Smith 
125416a05f60SBarry Smith    Level: beginner
125516a05f60SBarry Smith 
12561cc06b55SBarry Smith .seealso: [](ch_snes), `SNESFAS`, `PCMGSetType()`, `PCMGType`
1257a055c8caSBarry Smith E*/
12589371c9d4SSatish Balay typedef enum {
12599371c9d4SSatish Balay   SNES_FAS_MULTIPLICATIVE,
12609371c9d4SSatish Balay   SNES_FAS_ADDITIVE,
12619371c9d4SSatish Balay   SNES_FAS_FULL,
12629371c9d4SSatish Balay   SNES_FAS_KASKADE
12639371c9d4SSatish Balay } SNESFASType;
1264a055c8caSBarry Smith PETSC_EXTERN const char *const SNESFASTypes[];
1265a055c8caSBarry Smith 
1266a055c8caSBarry Smith /* called on the finest level FAS instance*/
1267a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetType(SNES, SNESFASType);
1268a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetType(SNES, SNESFASType *);
1269a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLevels(SNES, PetscInt, MPI_Comm *);
1270a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetLevels(SNES, PetscInt *);
1271a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCycleSNES(SNES, PetscInt, SNES *);
1272a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothUp(SNES, PetscInt);
1273a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothDown(SNES, PetscInt);
1274a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetCycles(SNES, PetscInt);
1275d142ab34SLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESFASSetMonitor(SNES, PetscViewerAndFormat *, PetscBool);
1276a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLog(SNES, PetscBool);
1277a055c8caSBarry Smith 
1278a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetGalerkin(SNES, PetscBool);
1279a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetGalerkin(SNES, PetscBool *);
128025acbd8eSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESFASGalerkinFunctionDefault(SNES, Vec, Vec, void *);
1281a055c8caSBarry Smith 
1282a055c8caSBarry Smith /* called on any level -- "Cycle" FAS instance */
1283a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmoother(SNES, SNES *);
1284a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherUp(SNES, SNES *);
1285a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherDown(SNES, SNES *);
1286a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetCorrection(SNES, SNES *);
1287a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInterpolation(SNES, Mat *);
1288a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRestriction(SNES, Mat *);
1289a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInjection(SNES, Mat *);
1290a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRScale(SNES, Vec *);
1291a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleSetCycles(SNES, PetscInt);
1292a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleIsFine(SNES, PetscBool *);
1293a055c8caSBarry Smith 
1294a055c8caSBarry Smith /* called on the (outer) finest level FAS to set/get parameters on any level instance */
1295a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInterpolation(SNES, PetscInt, Mat);
1296a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInterpolation(SNES, PetscInt, Mat *);
1297a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRestriction(SNES, PetscInt, Mat);
1298a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRestriction(SNES, PetscInt, Mat *);
1299a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInjection(SNES, PetscInt, Mat);
1300a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInjection(SNES, PetscInt, Mat *);
1301a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRScale(SNES, PetscInt, Vec);
1302a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRScale(SNES, PetscInt, Vec *);
1303020631bcSPeter Brune PETSC_EXTERN PetscErrorCode SNESFASSetContinuation(SNES, PetscBool);
1304a055c8caSBarry Smith 
1305a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmoother(SNES, PetscInt, SNES *);
1306a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherUp(SNES, PetscInt, SNES *);
1307a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherDown(SNES, PetscInt, SNES *);
1308a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCoarseSolve(SNES, SNES *);
1309a055c8caSBarry Smith 
1310a055c8caSBarry Smith /* parameters for full FAS */
1311a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASFullSetDownSweep(SNES, PetscBool);
1312a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCreateCoarseVec(SNES, Vec *);
1313a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASRestrict(SNES, Vec, Vec);
13146dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullSetTotal(SNES, PetscBool);
13156dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullGetTotal(SNES, PetscBool *);
1316a055c8caSBarry Smith 
13172eace19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetSNESVariableBounds(DM, SNES);
1318f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckDiscretization(SNES, DM, PetscReal, Vec, PetscReal, PetscReal[]);
1319f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckResidual(SNES, DM, Vec, PetscReal, PetscReal *);
1320f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckJacobian(SNES, DM, Vec, PetscReal, PetscBool *, PetscReal *);
13217f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckFromOptions(SNES, Vec);
13228e3a2eefSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESComputeJacobianAction(DM, Vec, Vec, Vec, void *);
13238e3a2eefSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCreateJacobianMF(DM, Vec, void *, Mat *);
132497276fddSZach Atkins 
132597276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALSetFunction(SNES, SNESFunctionFn *, void *ctx);
132697276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALGetFunction(SNES, SNESFunctionFn **, void **ctx);
132797276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALComputeFunction(SNES, Vec, Vec);
132897276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALGetLoadParameter(SNES, PetscReal *);
132997276fddSZach Atkins 
133097276fddSZach Atkins /*MC
133197276fddSZach Atkins    SNESNewtonALCorrectionType - the approach used by `SNESNEWTONAL` to determine
133297276fddSZach Atkins    the correction to the current increment. While the exact correction satisfies
133397276fddSZach Atkins    the constraint surface at every iteration, it also requires solving a quadratic
133497276fddSZach Atkins    equation which may not have real roots. Conversely, the normal correction is more
133597276fddSZach Atkins    efficient and always yields a real correction and is the default.
133697276fddSZach Atkins 
133797276fddSZach Atkins    Values:
133897276fddSZach Atkins +   `SNES_NEWTONAL_CORRECTION_EXACT` - choose the correction which exactly satisfies the constraint
1339d7c1f440SPierre Jolivet -   `SNES_NEWTONAL_CORRECTION_NORMAL` - choose the correction in the updated normal hyper-surface to the constraint surface
134097276fddSZach Atkins 
134197276fddSZach Atkins    Options Database Key:
134297276fddSZach Atkins . -snes_newtonal_correction_type <exact> - select type from <exact,normal>
134397276fddSZach Atkins 
134497276fddSZach Atkins    Level: intermediate
134597276fddSZach Atkins 
134697276fddSZach Atkins .seealso: `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetCorrectionType()`
134797276fddSZach Atkins M*/
134897276fddSZach Atkins typedef enum {
134997276fddSZach Atkins   SNES_NEWTONAL_CORRECTION_EXACT  = 0,
135097276fddSZach Atkins   SNES_NEWTONAL_CORRECTION_NORMAL = 1,
135197276fddSZach Atkins } SNESNewtonALCorrectionType;
135297276fddSZach Atkins PETSC_EXTERN const char *const SNESNewtonALCorrectionTypes[];
135397276fddSZach Atkins 
135497276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALSetCorrectionType(SNES, SNESNewtonALCorrectionType);
1355