xref: /petsc/include/petscsnes.h (revision d5def619400c53ceafc7498a3d382ce0d1e29f7d)
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 
2389804daf3SBarry Smith #include <petscdrawtypes.h>
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMonitorLGRange(SNES, PetscInt, PetscReal, void *);
240eafb4bcbSBarry Smith 
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetApplicationContext(SNES, void *);
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetApplicationContext(SNES, void *);
24349abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode SNESSetComputeApplicationContext(SNES, PetscErrorCode (*)(SNES, void **), PetscCtxDestroyFn *);
244184914b5SBarry Smith 
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESPythonSetType(SNES, const char[]);
246ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode SNESPythonGetType(SNES, const char *[]);
2471d6018f0SLisandro Dalcin 
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetFunctionDomainError(SNES);
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetFunctionDomainError(SNES, PetscBool *);
250cc6b0f04SFande Kong PETSC_EXTERN PetscErrorCode SNESGetJacobianDomainError(SNES, PetscBool *);
251cc6b0f04SFande Kong PETSC_EXTERN PetscErrorCode SNESSetJacobianDomainError(SNES);
252b351a90bSFande Kong PETSC_EXTERN PetscErrorCode SNESSetCheckJacobianDomainError(SNES, PetscBool);
2538383d7d7SFande Kong PETSC_EXTERN PetscErrorCode SNESGetCheckJacobianDomainError(SNES, PetscBool *);
2546a388c36SPeter Brune 
255edd03b47SJacob Faibussowitsch #define SNES_CONVERGED_TR_DELTA_DEPRECATED SNES_CONVERGED_TR_DELTA PETSC_DEPRECATED_ENUM(3, 12, 0, "SNES_DIVERGED_TR_DELTA", )
256435da068SBarry Smith /*E
257dc4c0fb0SBarry Smith     SNESConvergedReason - reason a `SNESSolve()` was determined to have converged or diverged
258435da068SBarry Smith 
259dc4c0fb0SBarry Smith    Values:
260dc4c0fb0SBarry Smith +  `SNES_CONVERGED_FNORM_ABS`      - 2-norm(F) <= abstol
261dc4c0fb0SBarry Smith .  `SNES_CONVERGED_FNORM_RELATIVE` - 2-norm(F) <= rtol*2-norm(F(x_0)) where x_0 is the initial guess
262dc4c0fb0SBarry Smith .  `SNES_CONVERGED_SNORM_RELATIVE` - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
263379ef8d2SAlexander .  `SNES_CONVERGED_USER`           - The user has indicated convergence for an arbitrary reason
264dc4c0fb0SBarry Smith .  `SNES_DIVERGED_FUNCTION_COUNT`  - The user provided function has been called more times than the maximum set in `SNESSetTolerances()`
265dc4c0fb0SBarry Smith .  `SNES_DIVERGED_DTOL`            - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
266dc4c0fb0SBarry Smith .  `SNES_DIVERGED_FNORM_NAN`       - the 2-norm of the current function evaluation is not-a-number (NaN), this
267dc4c0fb0SBarry Smith                                      is usually caused by a division of 0 by 0.
268dc4c0fb0SBarry Smith .  `SNES_DIVERGED_MAX_IT`          - `SNESSolve()` has reached the maximum number of iterations requested
269dc4c0fb0SBarry Smith .  `SNES_DIVERGED_LINE_SEARCH`     - The line search has failed. This only occurs for `SNES` solvers that use a line search
270dc4c0fb0SBarry Smith .  `SNES_DIVERGED_LOCAL_MIN`       - the algorithm seems to have stagnated at a local minimum that is not zero.
271379ef8d2SAlexander .  `SNES_DIVERGED_USER`            - The user has indicated divergence for an arbitrary reason
272379ef8d2SAlexander -  `SNES_CONVERGED_ITERATING       - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
273f203c74bSBarry Smith 
27416a05f60SBarry Smith    Level: beginner
27516a05f60SBarry Smith 
276dc4c0fb0SBarry Smith     Notes:
2770b4b7b1cSBarry 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
2780b4b7b1cSBarry Smith    (in this case we recommend
279dc4c0fb0SBarry Smith    testing with `-pc_type lu` to eliminate the linear solver as the cause of the problem).
280dc4c0fb0SBarry Smith 
2810b4b7b1cSBarry Smith    `SNES_DIVERGED_LOCAL_MIN` can only occur when using a `SNES` solver that uses a line search (`SNESLineSearch`).
282f203c74bSBarry Smith    The line search wants to minimize Q(alpha) = 1/2 || F(x + alpha s) ||^2_2  this occurs
283f203c74bSBarry 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
284f203c74bSBarry Smith    you get Q'(alpha) = -F(x)^T F'(x)^(-1)^T F'(x+alpha s)F(x+alpha s); when alpha = 0
285f203c74bSBarry Smith    Q'(0) = - ||F(x)||^2_2 which is always NEGATIVE if F'(x) is invertible. This means the Newton
286f203c74bSBarry Smith    direction is a descent direction and the line search should succeed if alpha is small enough.
287f203c74bSBarry Smith 
288f203c74bSBarry Smith    If F'(x) is NOT invertible AND F'(x)^T F(x) = 0 then Q'(0) = 0 and the Newton direction
289f203c74bSBarry Smith    is NOT a descent direction so the line search will fail. All one can do at this point
290f203c74bSBarry Smith    is change the initial guess and try again.
291f203c74bSBarry Smith 
292f203c74bSBarry Smith    An alternative explanation: Newton's method can be regarded as replacing the function with
293f203c74bSBarry Smith    its linear approximation and minimizing the 2-norm of that. That is F(x+s) approx F(x) + F'(x)s
294f203c74bSBarry Smith    so we minimize || F(x) + F'(x) s ||^2_2; do this using Least Squares. If F'(x) is invertible then
295f203c74bSBarry 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
2966aad120cSJose E. Roman    exists a nontrivial (that is F'(x)s != 0) solution to the equation and this direction is
297f203c74bSBarry 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)
298f203c74bSBarry 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
299f203c74bSBarry Smith    and F'(x)^T F'(x) has no negative eigenvalues Q'(0) < 0 so s is a descent direction and the line
300f203c74bSBarry Smith    search should succeed for small enough alpha.
301f203c74bSBarry Smith 
302f203c74bSBarry Smith    Note that this RARELY happens in practice. Far more likely the linear system is not being solved
303f203c74bSBarry Smith    (well enough?) or the Jacobian is wrong.
304f203c74bSBarry Smith 
30587497f52SBarry Smith    `SNES_DIVERGED_MAX_IT` means that the solver reached the maximum number of iterations without satisfying any
30687497f52SBarry Smith    convergence criteria. `SNES_CONVERGED_ITS` means that `SNESConvergedSkip()` was chosen as the convergence test;
3074d0a8057SBarry Smith    thus the usual convergence criteria have not been checked and may or may not be satisfied.
308f203c74bSBarry Smith 
3091cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`, `SNESSetTolerances()`
310435da068SBarry Smith E*/
311184914b5SBarry Smith typedef enum {                       /* converged */
31201b82886SBarry Smith   SNES_CONVERGED_FNORM_ABS      = 2, /* ||F|| < atol */
31301b82886SBarry Smith   SNES_CONVERGED_FNORM_RELATIVE = 3, /* ||F|| < rtol*||F_initial|| */
3145358d0d4SBarry Smith   SNES_CONVERGED_SNORM_RELATIVE = 4, /* Newton computed step size small; || delta x || < stol || x ||*/
3153f149594SLisandro Dalcin   SNES_CONVERGED_ITS            = 5, /* maximum iterations reached */
316379ef8d2SAlexander   SNES_BREAKOUT_INNER_ITER      = 6, /* Flag to break out of inner loop after checking custom convergence.
317379ef8d2SAlexander                                         It is used in multi-phase flow when state changes */
318379ef8d2SAlexander   SNES_CONVERGED_USER           = 7, /* The user has indicated convergence for an arbitrary reason */
319184914b5SBarry Smith   /* diverged */
32046a9e3ceSBarry Smith   SNES_DIVERGED_FUNCTION_DOMAIN      = -1, /* the new x location passed the function is not in the domain of F */
321184914b5SBarry Smith   SNES_DIVERGED_FUNCTION_COUNT       = -2,
32246a9e3ceSBarry Smith   SNES_DIVERGED_LINEAR_SOLVE         = -3, /* the linear solve failed */
323184914b5SBarry Smith   SNES_DIVERGED_FNORM_NAN            = -4,
324184914b5SBarry Smith   SNES_DIVERGED_MAX_IT               = -5,
325647a2e1fSBarry Smith   SNES_DIVERGED_LINE_SEARCH          = -6,  /* the line search failed */
3261e633543SBarry Smith   SNES_DIVERGED_INNER                = -7,  /* inner solve failed */
32758c775ebSBarry Smith   SNES_DIVERGED_LOCAL_MIN            = -8,  /* || J^T b || is small, implies converged to local minimum of F() */
328e37c518bSBarry Smith   SNES_DIVERGED_DTOL                 = -9,  /* || F || > divtol*||F_initial|| */
32907b62357SFande Kong   SNES_DIVERGED_JACOBIAN_DOMAIN      = -10, /* Jacobian calculation does not make sense */
3301c6b2ff8SBarry Smith   SNES_DIVERGED_TR_DELTA             = -11,
331c78072a7SPatrick Sanan   SNES_CONVERGED_TR_DELTA_DEPRECATED = -11,
332379ef8d2SAlexander   SNES_DIVERGED_USER                 = -12, /* The user has indicated divergence for an arbitrary reason */
333e37c518bSBarry Smith 
3349371c9d4SSatish Balay   SNES_CONVERGED_ITERATING = 0
3359371c9d4SSatish Balay } SNESConvergedReason;
336014dd563SJed Brown PETSC_EXTERN const char *const *SNESConvergedReasons;
337184914b5SBarry Smith 
338c838673bSBarry Smith /*MC
339c838673bSBarry Smith    SNES_CONVERGED_FNORM_ABS - 2-norm(F) <= abstol
340c838673bSBarry Smith 
341c838673bSBarry Smith    Level: beginner
342c838673bSBarry Smith 
3431cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
344c838673bSBarry Smith M*/
345c838673bSBarry Smith 
346c838673bSBarry Smith /*MC
347c838673bSBarry Smith    SNES_CONVERGED_FNORM_RELATIVE - 2-norm(F) <= rtol*2-norm(F(x_0)) where x_0 is the initial guess
348c838673bSBarry Smith 
349c838673bSBarry Smith    Level: beginner
350c838673bSBarry Smith 
3511cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
352c838673bSBarry Smith M*/
353c838673bSBarry Smith 
354c838673bSBarry Smith /*MC
355c60f73f4SPeter Brune   SNES_CONVERGED_SNORM_RELATIVE - The 2-norm of the last step <= stol * 2-norm(x) where x is the current
35687497f52SBarry Smith   solution and stol is the 4th argument to `SNESSetTolerances()`
357c838673bSBarry Smith 
358af27ebaaSBarry Smith   Options Database Key:
359ae9be289SBarry Smith   -snes_stol <stol> - the step tolerance
360ae9be289SBarry Smith 
361c838673bSBarry Smith    Level: beginner
362c838673bSBarry Smith 
3631cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
364c838673bSBarry Smith M*/
365c838673bSBarry Smith 
366c838673bSBarry Smith /*MC
367c838673bSBarry Smith    SNES_DIVERGED_FUNCTION_COUNT - The user provided function has been called more times then the final
36887497f52SBarry Smith    argument to `SNESSetTolerances()`
369c838673bSBarry Smith 
370c838673bSBarry Smith    Level: beginner
371c838673bSBarry Smith 
3721cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
373c838673bSBarry Smith M*/
374c838673bSBarry Smith 
375c838673bSBarry Smith /*MC
37687497f52SBarry Smith    SNES_DIVERGED_DTOL - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
377e37c518bSBarry Smith 
378e37c518bSBarry Smith    Level: beginner
379e37c518bSBarry Smith 
3801cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESSetDivergenceTolerance()`
381e37c518bSBarry Smith M*/
382e37c518bSBarry Smith 
383e37c518bSBarry Smith /*MC
384c838673bSBarry Smith    SNES_DIVERGED_FNORM_NAN - the 2-norm of the current function evaluation is not-a-number (NaN), this
385c838673bSBarry Smith    is usually caused by a division of 0 by 0.
386c838673bSBarry Smith 
387c838673bSBarry Smith    Level: beginner
388c838673bSBarry Smith 
3891cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
390c838673bSBarry Smith M*/
391c838673bSBarry Smith 
392c838673bSBarry Smith /*MC
393c838673bSBarry Smith    SNES_DIVERGED_MAX_IT - SNESSolve() has reached the maximum number of iterations requested
394c838673bSBarry Smith 
395c838673bSBarry Smith    Level: beginner
396c838673bSBarry Smith 
3971cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
398c838673bSBarry Smith M*/
399c838673bSBarry Smith 
400c838673bSBarry Smith /*MC
40187497f52SBarry Smith    SNES_DIVERGED_LINE_SEARCH - The line search has failed. This only occurs for a `SNES` solvers that use a line search
402c838673bSBarry Smith 
403c838673bSBarry Smith    Level: beginner
404c838673bSBarry Smith 
4051cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESLineSearch`
406c838673bSBarry Smith M*/
407c838673bSBarry Smith 
408c838673bSBarry Smith /*MC
40946a9e3ceSBarry Smith    SNES_DIVERGED_LOCAL_MIN - the algorithm seems to have stagnated at a local minimum that is not zero.
41087497f52SBarry Smith    See the manual page for `SNESConvergedReason` for more details
411c838673bSBarry Smith 
412c838673bSBarry Smith    Level: beginner
413c838673bSBarry Smith 
4141cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
415c838673bSBarry Smith M*/
416c838673bSBarry Smith 
417c838673bSBarry Smith /*MC
41887497f52SBarry Smith    SNES_CONERGED_ITERATING - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
419c838673bSBarry Smith 
420c838673bSBarry Smith    Level: beginner
421c838673bSBarry Smith 
4221cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
423c838673bSBarry Smith M*/
424c838673bSBarry Smith 
425014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetConvergenceTest(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *), void *, PetscErrorCode (*)(void *));
4268d359177SBarry Smith PETSC_EXTERN PetscErrorCode SNESConvergedDefault(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
427e2a6519dSDmitry Karpeev PETSC_EXTERN PetscErrorCode SNESConvergedSkip(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
428649ef022SMatthew Knepley PETSC_EXTERN PetscErrorCode SNESConvergedCorrectPressure(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
429014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetConvergedReason(SNES, SNESConvergedReason *);
430c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode SNESGetConvergedReasonString(SNES, const char **);
43133866048SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetConvergedReason(SNES, SNESConvergedReason);
432ddbbdb52SLois Curfman McInnes 
433edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "SNESConvergedSkip()", ) static inline void SNESSkipConverged(void)
434d71ae5a4SJacob Faibussowitsch { /* never called */
4359371c9d4SSatish Balay }
4368ea1b3e6SJed Brown #define SNESSkipConverged (SNESSkipConverged, SNESConvergedSkip)
43711f088b5SMatthew G Knepley 
4389bcc50f1SBarry Smith /*S
4398434afd1SBarry Smith   SNESInitialGuessFn - A prototype of a `SNES` compute initial guess function that would be passed to `SNESSetComputeInitialGuess()`
4409bcc50f1SBarry Smith 
4419bcc50f1SBarry Smith   Calling Sequence:
4429bcc50f1SBarry Smith + snes  - `SNES` context
4439bcc50f1SBarry Smith . u   - output vector to contain initial guess
4449bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4459bcc50f1SBarry Smith 
4469bcc50f1SBarry Smith   Level: beginner
4479bcc50f1SBarry Smith 
4488434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetComputeInitialGuess()`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESFunctionFn`
4499bcc50f1SBarry Smith S*/
4508434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESInitialGuessFn)(SNES snes, Vec u, void *ctx);
4519bcc50f1SBarry Smith 
4529bcc50f1SBarry Smith /*S
4538434afd1SBarry Smith   SNESFunctionFn - A prototype of a `SNES` evaluation function that would be passed to `SNESSetFunction()`
4549bcc50f1SBarry Smith 
4559bcc50f1SBarry Smith   Calling Sequence:
4569bcc50f1SBarry Smith + snes  - `SNES` context
4579bcc50f1SBarry Smith . u   - input vector
4589bcc50f1SBarry Smith . F   - function vector
4599bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4609bcc50f1SBarry Smith 
4619bcc50f1SBarry Smith   Level: beginner
4629bcc50f1SBarry Smith 
4638434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
4649bcc50f1SBarry Smith S*/
4658434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESFunctionFn)(SNES snes, Vec u, Vec F, void *ctx);
4669bcc50f1SBarry Smith 
4679bcc50f1SBarry Smith /*S
4688434afd1SBarry Smith   SNESObjectiveFn - A prototype of a `SNES` objective evaluation function that would be passed to `SNESSetObjective()`
4699bcc50f1SBarry Smith 
4709bcc50f1SBarry Smith   Calling Sequence:
4719bcc50f1SBarry Smith + snes  - `SNES` context
4729bcc50f1SBarry Smith . u   - input vector
4739bcc50f1SBarry Smith . o   - output value
4749bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4759bcc50f1SBarry Smith 
4769bcc50f1SBarry Smith   Level: beginner
4779bcc50f1SBarry Smith 
4788434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
4799bcc50f1SBarry Smith S*/
4808434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESObjectiveFn)(SNES snes, Vec u, PetscReal *o, void *ctx);
4819bcc50f1SBarry Smith 
4829bcc50f1SBarry Smith /*S
4838434afd1SBarry Smith   SNESJacobianFn - A prototype of a `SNES` Jacobian evaluation function that would be passed to `SNESSetJacobian()`
4849bcc50f1SBarry Smith 
4859bcc50f1SBarry Smith   Calling Sequence:
4869bcc50f1SBarry Smith + snes   - the `SNES` context obtained from `SNESCreate()`
4879bcc50f1SBarry Smith . u    - input vector
4889bcc50f1SBarry Smith . Amat - (approximate) Jacobian matrix
4899bcc50f1SBarry Smith . Pmat - matrix used to construct the preconditioner, often the same as `Amat`
4909bcc50f1SBarry Smith - ctx  - [optional] user-defined context for matrix evaluation routine
4919bcc50f1SBarry Smith 
4929bcc50f1SBarry Smith   Level: beginner
4939bcc50f1SBarry Smith 
4948434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESNGSFn`
4959bcc50f1SBarry Smith S*/
4968434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESJacobianFn)(SNES snes, Vec u, Mat Amat, Mat Pmat, void *ctx);
4979bcc50f1SBarry Smith 
4989bcc50f1SBarry Smith /*S
4998434afd1SBarry Smith   SNESNGSFn - A prototype of a `SNES` nonlinear Gauss-Seidel function that would be passed to `SNESSetNGS()`
5009bcc50f1SBarry Smith 
5019bcc50f1SBarry Smith   Calling Sequence:
5029bcc50f1SBarry Smith + snes   - the `SNES` context obtained from `SNESCreate()`
5039bcc50f1SBarry Smith . u    - the current solution, updated in place
504dd8e379bSPierre Jolivet . b    - the right-hand side vector (which may be `NULL`)
5059bcc50f1SBarry Smith - ctx  - [optional] user-defined context for matrix evaluation routine
5069bcc50f1SBarry Smith 
5079bcc50f1SBarry Smith   Level: beginner
5089bcc50f1SBarry Smith 
5098434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`
5109bcc50f1SBarry Smith S*/
5118434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESNGSFn)(SNES snes, Vec u, Vec b, void *ctx);
5129bcc50f1SBarry Smith 
51353e5d35bSStefano Zampini /*S
51453e5d35bSStefano Zampini   SNESUpdateFn - A prototype of a `SNES` update function that would be passed to `SNESSetUpdate()`
51553e5d35bSStefano Zampini 
51653e5d35bSStefano Zampini   Calling Sequence:
51753e5d35bSStefano Zampini + snes - `SNES` context
51853e5d35bSStefano Zampini - step - the current iteration index
51953e5d35bSStefano Zampini 
52053e5d35bSStefano Zampini   Level: advanced
52153e5d35bSStefano Zampini 
52253e5d35bSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESSetUpdate()`
52353e5d35bSStefano Zampini S*/
52453e5d35bSStefano Zampini PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESUpdateFn)(SNES snes, PetscInt step);
52553e5d35bSStefano Zampini 
526b67197daSBarry Smith /* --------- Solving systems of nonlinear equations --------------- */
5278434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetFunction(SNES, Vec, SNESFunctionFn *, void *);
5288434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetFunction(SNES, Vec *, SNESFunctionFn **, void **);
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESComputeFunction(SNES, Vec, Vec);
530bbc1464cSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeMFFunction(SNES, Vec, Vec);
53125acbd8eSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESSetInitialFunction(SNES, Vec);
53232f3f7c2SPeter Brune 
5338434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetJacobian(SNES, Mat, Mat, SNESJacobianFn *, void *);
5348434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetJacobian(SNES, Mat *, Mat *, SNESJacobianFn **, void **);
5358434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESObjectiveComputeFunctionDefaultFD;
5368434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefault;
5378434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefaultColor;
53878e7fe0eSHong Zhang PETSC_EXTERN PetscErrorCode SNESPruneJacobianColor(SNES, Mat, Mat);
5398434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetComputeInitialGuess(SNES, SNESInitialGuessFn *, void *);
5408434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetPicard(SNES, Vec, SNESFunctionFn *, Mat, Mat, SNESJacobianFn *, void *);
5418434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetPicard(SNES, Vec *, SNESFunctionFn **, Mat *, Mat *, SNESJacobianFn **, void **);
5428434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeFunction;
5438434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeMFFunction;
5448434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESPicardComputeJacobian;
54517bae607SBarry Smith 
5468434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetObjective(SNES, SNESObjectiveFn *, void *);
5478434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetObjective(SNES, SNESObjectiveFn **, void **);
5482a4ee8f2SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeObjective(SNES, Vec, PetscReal *);
5492a4ee8f2SPeter Brune 
55053e5d35bSStefano Zampini PETSC_EXTERN PetscErrorCode SNESSetUpdate(SNES, SNESUpdateFn *);
55153e5d35bSStefano Zampini 
552534ebe21SPeter Brune /*E
55316a05f60SBarry Smith    SNESNormSchedule - Frequency with which the norm is computed during a nonliner solve
554534ebe21SPeter Brune 
555dc4c0fb0SBarry Smith    Values:
556dc4c0fb0SBarry Smith +   `SNES_NORM_DEFAULT`            - use the default behavior for the current `SNESType`
557dc4c0fb0SBarry Smith .   `SNES_NORM_NONE`               - avoid all norm computations
558dc4c0fb0SBarry Smith .   `SNES_NORM_ALWAYS`             - compute the norms whenever possible
559dc4c0fb0SBarry Smith .   `SNES_NORM_INITIAL_ONLY`       - compute the norm only when the algorithm starts
560dc4c0fb0SBarry Smith .   `SNES_NORM_FINAL_ONLY`         - compute the norm only when the algorithm finishes
561dc4c0fb0SBarry Smith -   `SNES_NORM_INITIAL_FINAL_ONLY` - compute the norm at the start and end of the algorithm
562534ebe21SPeter Brune 
56316a05f60SBarry Smith    Level: advanced
56416a05f60SBarry Smith 
565534ebe21SPeter Brune    Notes:
566dc4c0fb0SBarry Smith    Support for these is highly dependent on the solver.
567dc4c0fb0SBarry Smith 
568dc4c0fb0SBarry Smith    Some options limit the convergence tests that can be used.
569dc4c0fb0SBarry Smith 
570dc4c0fb0SBarry Smith    The `SNES_NORM_NONE` option is most commonly used when the nonlinear solver is being used as a smoother, for example for `SNESFAS`
571dc4c0fb0SBarry Smith 
572534ebe21SPeter Brune    This is primarily used to turn off extra norm and function computation
573534ebe21SPeter Brune    when the solvers are composed.
574534ebe21SPeter Brune 
5751cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
576db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
577534ebe21SPeter Brune E*/
5789371c9d4SSatish Balay typedef enum {
5799371c9d4SSatish Balay   SNES_NORM_DEFAULT            = -1,
580fdacfa88SPeter Brune   SNES_NORM_NONE               = 0,
581365a6726SPeter Brune   SNES_NORM_ALWAYS             = 1,
582fdacfa88SPeter Brune   SNES_NORM_INITIAL_ONLY       = 2,
583fdacfa88SPeter Brune   SNES_NORM_FINAL_ONLY         = 3,
5849371c9d4SSatish Balay   SNES_NORM_INITIAL_FINAL_ONLY = 4
5859371c9d4SSatish Balay } SNESNormSchedule;
586365a6726SPeter Brune PETSC_EXTERN const char *const *const SNESNormSchedules;
5871957e957SBarry Smith 
588534ebe21SPeter Brune /*MC
58916a05f60SBarry Smith    SNES_NORM_NONE - Don't compute function and its L2 norm when possible
590534ebe21SPeter Brune 
591534ebe21SPeter Brune    Level: advanced
592534ebe21SPeter Brune 
593dc4c0fb0SBarry Smith    Note:
594534ebe21SPeter Brune    This is most useful for stationary solvers with a fixed number of iterations used as smoothers.
595534ebe21SPeter Brune 
5961cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_DEFAULT`
597534ebe21SPeter Brune M*/
598534ebe21SPeter Brune 
599fdacfa88SPeter Brune /*MC
600365a6726SPeter Brune    SNES_NORM_ALWAYS - Compute the function and its L2 norm at each iteration.
601534ebe21SPeter Brune 
602fdacfa88SPeter Brune    Level: advanced
603fdacfa88SPeter Brune 
604dc4c0fb0SBarry Smith    Note:
605fdacfa88SPeter Brune    Most solvers will use this no matter what norm type is passed to them.
606fdacfa88SPeter Brune 
6071cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_NONE`
608fdacfa88SPeter Brune M*/
609534ebe21SPeter Brune 
610534ebe21SPeter Brune /*MC
611534ebe21SPeter Brune    SNES_NORM_INITIAL_ONLY - Compute the function and its L2 at iteration 0, but do not update it.
612534ebe21SPeter Brune 
613534ebe21SPeter Brune    Level: advanced
614534ebe21SPeter Brune 
615534ebe21SPeter Brune    Notes:
616dc4c0fb0SBarry Smith    This method is useful in composed methods, when a true solution might actually be found before `SNESSolve()` is called.
617534ebe21SPeter Brune    This option enables the solve to abort on the zeroth iteration if this is the case.
618534ebe21SPeter Brune 
619534ebe21SPeter Brune    For solvers that require the computation of the L2 norm of the function as part of the method, this merely cancels
620534ebe21SPeter Brune    the norm computation at the last iteration (if possible).
621534ebe21SPeter Brune 
6221cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_FINAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
623534ebe21SPeter Brune M*/
624534ebe21SPeter Brune 
625534ebe21SPeter Brune /*MC
626534ebe21SPeter Brune    SNES_NORM_FINAL_ONLY - Compute the function and its L2 norm on only the final iteration.
627534ebe21SPeter Brune 
628534ebe21SPeter Brune    Level: advanced
629534ebe21SPeter Brune 
630dc4c0fb0SBarry Smith    Note:
631534ebe21SPeter Brune    For solvers that require the computation of the L2 norm of the function as part of the method, behaves
63287497f52SBarry Smith    exactly as `SNES_NORM_DEFAULT`.  This method is useful when the function is gotten after `SNESSolve()` and
633534ebe21SPeter Brune    used in subsequent computation for methods that do not need the norm computed during the rest of the
634534ebe21SPeter Brune    solution procedure.
635534ebe21SPeter Brune 
6361cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_INITIAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
637534ebe21SPeter Brune M*/
638534ebe21SPeter Brune 
639534ebe21SPeter Brune /*MC
640534ebe21SPeter Brune    SNES_NORM_INITIAL_FINAL_ONLY - Compute the function and its L2 norm on only the initial and final iterations.
641534ebe21SPeter Brune 
642534ebe21SPeter Brune    Level: advanced
643534ebe21SPeter Brune 
644dc4c0fb0SBarry Smith    Note:
64587497f52SBarry Smith    This method combines the benefits of `SNES_NORM_INITIAL_ONLY` and `SNES_NORM_FINAL_ONLY`.
646534ebe21SPeter Brune 
6471cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_SNES_NORM_INITIAL_ONLY`, `SNES_NORM_FINAL_ONLY`
648534ebe21SPeter Brune M*/
649534ebe21SPeter Brune 
650365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetNormSchedule(SNES, SNESNormSchedule);
651365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetNormSchedule(SNES, SNESNormSchedule *);
652c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetFunctionNorm(SNES, PetscReal);
653c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESGetFunctionNorm(SNES, PetscReal *);
654c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetUpdateNorm(SNES, PetscReal *);
655c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetSolutionNorm(SNES, PetscReal *);
656534ebe21SPeter Brune 
65747073ea2SPeter Brune /*E
65847073ea2SPeter Brune    SNESFunctionType - Type of function computed
65947073ea2SPeter Brune 
660dc4c0fb0SBarry Smith    Values:
661dc4c0fb0SBarry Smith +  `SNES_FUNCTION_DEFAULT`          - the default behavior for the current `SNESType`
662dc4c0fb0SBarry Smith .  `SNES_FUNCTION_UNPRECONDITIONED` - the original function provided
663dc4c0fb0SBarry Smith -  `SNES_FUNCTION_PRECONDITIONED`   - the modification of the function by the preconditioner
66447073ea2SPeter Brune 
66516a05f60SBarry Smith    Level: advanced
66616a05f60SBarry Smith 
667dc4c0fb0SBarry Smith    Note:
668dc4c0fb0SBarry Smith    Support for these is dependent on the solver.
669dc4c0fb0SBarry Smith 
6701cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
671db781477SPatrick Sanan           `KSPSetConvergenceTest()`, `KSPSetPCSide()`
67247073ea2SPeter Brune E*/
6739371c9d4SSatish Balay typedef enum {
6749371c9d4SSatish Balay   SNES_FUNCTION_DEFAULT          = -1,
67547073ea2SPeter Brune   SNES_FUNCTION_UNPRECONDITIONED = 0,
6769371c9d4SSatish Balay   SNES_FUNCTION_PRECONDITIONED   = 1
6779371c9d4SSatish Balay } SNESFunctionType;
67847073ea2SPeter Brune PETSC_EXTERN const char *const *const SNESFunctionTypes;
67947073ea2SPeter Brune 
68047073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetFunctionType(SNES, SNESFunctionType);
68147073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetFunctionType(SNES, SNESFunctionType *);
682f1c6b773SPeter Brune 
6838434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNGS(SNES, SNESNGSFn *, void *);
6848434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNGS(SNES, SNESNGSFn **, void **);
685be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGS(SNES, Vec, Vec);
686b6266c6eSPeter Brune 
687be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetSweeps(SNES, PetscInt);
688be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetSweeps(SNES, PetscInt *);
689be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt);
690be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
691badc63e7SPeter Brune 
6924fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES, PetscBool);
6934fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES, PetscBool *);
6944fc747eaSLawrence Mitchell 
6953ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode SNESShellGetContext(SNES, void *);
696014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESShellSetContext(SNES, void *);
697014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESShellSetSolve(SNES, PetscErrorCode (*)(SNES, Vec));
698646217ecSPeter Brune 
699c5ae4b9aSBarry Smith /* --------- Routines specifically for line search methods --------------- */
700c5ae4b9aSBarry Smith 
701872b6db9SPeter Brune /*S
70287497f52SBarry Smith    SNESLineSearch - Abstract PETSc object that manages line-search operations for nonlinear solvers
7039e764e56SPeter Brune 
7049e764e56SPeter Brune    Level: beginner
7059e764e56SPeter Brune 
7061cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNES`
7079e764e56SPeter Brune S*/
708907376e6SBarry Smith typedef struct _p_LineSearch *SNESLineSearch;
7099e764e56SPeter Brune 
7109e764e56SPeter Brune /*J
7110b4b7b1cSBarry Smith    SNESLineSearchType - String with the name of a PETSc line search method `SNESLineSearch`. Provides all the linesearches for the nonlinear solvers, `SNES`,
7120b4b7b1cSBarry Smith                         in PETSc.
7139e764e56SPeter Brune 
714ceaaa498SBarry Smith    Values:
715ceaaa498SBarry Smith +  `SNESLINESEARCHBASIC`   - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
716ceaaa498SBarry Smith .  `SNESLINESEARCHBT`      - Backtracking line search over the L2 norm of the function
717ceaaa498SBarry Smith .  `SNESLINESEARCHL2`      - Secant line search over the L2 norm of the function
7180b4b7b1cSBarry Smith .  `SNESLINESEARCHCP`      - Critical point secant line search assuming $F(x) = \nabla G(x)$ for some unknown $G(x)$
719ceaaa498SBarry Smith .  `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch
720ceaaa498SBarry Smith -  `SNESLINESEARCHSHELL`   - User provided `SNESLineSearch` implementation
721ceaaa498SBarry Smith 
72295bd0b28SBarry Smith    Level: beginner
72395bd0b28SBarry Smith 
7240b4b7b1cSBarry Smith    Note:
7250b4b7b1cSBarry Smith    Use `SNESLineSearchSetType()` or the options database key `-snes_linesearch_type` to set
7260b4b7b1cSBarry Smith    the specific line search algorithm to use with a given `SNES` object. Not all `SNESType` can utilize a line search.
7270b4b7b1cSBarry Smith 
7281cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetType()`, `SNES`
7299e764e56SPeter Brune J*/
73019fd82e9SBarry Smith typedef const char *SNESLineSearchType;
7319e764e56SPeter Brune #define SNESLINESEARCHBT        "bt"
732d4c6564cSPatrick Farrell #define SNESLINESEARCHNLEQERR   "nleqerr"
733c87759e9SPeter Brune #define SNESLINESEARCHBASIC     "basic"
7340b00b554SBarry Smith #define SNESLINESEARCHNONE      "none"
7359e764e56SPeter Brune #define SNESLINESEARCHL2        "l2"
7369e764e56SPeter Brune #define SNESLINESEARCHCP        "cp"
7379e764e56SPeter Brune #define SNESLINESEARCHSHELL     "shell"
738b5badacbSBarry Smith #define SNESLINESEARCHNCGLINEAR "ncglinear"
739b9d635d7SJonas Heinzmann #define SNESLINESEARCHBISECTION "bisection"
7409e764e56SPeter Brune 
741140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESList;
742014dd563SJed Brown PETSC_EXTERN PetscClassId      SNESLINESEARCH_CLASSID;
743140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESLineSearchList;
7449e764e56SPeter Brune 
745b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_LINEAR    1
746b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_QUADRATIC 2
747b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_CUBIC     3
7489e764e56SPeter Brune 
7499bcc50f1SBarry Smith /*S
7508434afd1SBarry Smith   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that projects a vector onto the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
7519bcc50f1SBarry Smith 
7529bcc50f1SBarry Smith   Calling Sequence:
7539bcc50f1SBarry Smith + snes  - `SNES` context
7549bcc50f1SBarry Smith - u     - the vector to project to the bounds
7559bcc50f1SBarry Smith 
7569bcc50f1SBarry Smith   Level: advanced
7579bcc50f1SBarry Smith 
7589bcc50f1SBarry Smith   Note:
7598434afd1SBarry Smith   The deprecated `SNESLineSearchVIProjectFunc` still works as a replacement for `SNESLineSearchVIProjectFn` *.
7609bcc50f1SBarry Smith 
7619bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
7629bcc50f1SBarry Smith S*/
7638434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESLineSearchVIProjectFn)(SNES snes, Vec u);
7648434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVIProjectFn *SNESLineSearchVIProjectFunc; // deprecated
7659bcc50f1SBarry Smith 
7669bcc50f1SBarry Smith /*S
7678434afd1SBarry Smith   SNESLineSearchVIProjectFn - A prototype of a `SNES` function that computes the norm of the active set variables in a vector in a VI solve,
7689bcc50f1SBarry Smith   passed to `SNESLineSearchSetVIFunctions()`
7699bcc50f1SBarry Smith 
7709bcc50f1SBarry Smith   Calling Sequence:
7719bcc50f1SBarry Smith + snes  - `SNES` context
7729bcc50f1SBarry Smith . f     - the vector to compute the norm of
7739bcc50f1SBarry Smith . u     - the current solution, entries that are on the VI bounds are ignored
7749bcc50f1SBarry Smith - fnorm - the resulting norm
7759bcc50f1SBarry Smith 
7769bcc50f1SBarry Smith   Level: advanced
7779bcc50f1SBarry Smith 
7789bcc50f1SBarry Smith   Note:
7798434afd1SBarry Smith   The deprecated `SNESLineSearchVINormFunc` still works as a replacement for `SNESLineSearchVINormFn` *.
7809bcc50f1SBarry Smith 
7819bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
7829bcc50f1SBarry Smith S*/
7838434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESLineSearchVINormFn)(SNES snes, Vec f, Vec u, PetscReal *fnorm);
7848434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVINormFn *SNESLineSearchVINormFunc; // deprecated
7859bcc50f1SBarry Smith 
786*d5def619SJonas Heinzmann /*S
787*d5def619SJonas Heinzmann   SNESLineSearchVIDirDerivFn - A prototype of a `SNES` function that computes the directional derivative considering the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
788*d5def619SJonas Heinzmann 
789*d5def619SJonas Heinzmann   Calling Sequence:
790*d5def619SJonas Heinzmann + snes  - `SNES` context
791*d5def619SJonas Heinzmann . f     - the function vector to compute the directional derivative with
792*d5def619SJonas Heinzmann . u     - the current solution, entries that are on the VI bounds are ignored
793*d5def619SJonas Heinzmann . y     - the direction to compute the directional derivative
794*d5def619SJonas Heinzmann - fty   - the resulting directional derivative
795*d5def619SJonas Heinzmann 
796*d5def619SJonas Heinzmann   Level: advanced
797*d5def619SJonas Heinzmann 
798*d5def619SJonas Heinzmann .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchVIProjectFn`, `SNESLineSearchVIProjectFn`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetVIFunctions()`
799*d5def619SJonas Heinzmann S*/
800*d5def619SJonas Heinzmann PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESLineSearchVIDirDerivFn)(SNES snes, Vec f, Vec u, Vec y, PetscScalar *fty);
801*d5def619SJonas Heinzmann 
8028434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESLineSearchApplyFn)(SNESLineSearch);
8038434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchApplyFn *SNESLineSearchApplyFunc; // deprecated
8048434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(SNESLineSearchShellApplyFn)(SNESLineSearch, void *);
8058434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchShellApplyFn *SNESLineSearchUserFunc; // deprecated
8069e764e56SPeter Brune 
807014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate(MPI_Comm, SNESLineSearch *);
808014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchReset(SNESLineSearch);
809014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchView(SNESLineSearch, PetscViewer);
810014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *);
811a80ff896SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetType(SNESLineSearch, SNESLineSearchType *);
81219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetType(SNESLineSearch, SNESLineSearchType);
813014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch);
814ed07d7d7SPeter Brune PETSC_EXTERN PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch, PetscErrorCode (*)(SNES, Vec, Vec));
815014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetUp(SNESLineSearch);
816014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchApply(SNESLineSearch, Vec, Vec, PetscReal *, Vec);
817014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch, Vec, Vec, PetscBool *);
818014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *);
819fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch, PetscInt);
8209e764e56SPeter Brune 
8219bd66eb0SPeter Brune /* set the functions for precheck and postcheck */
82286d74e61SPeter Brune 
823f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx);
824f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx);
8259e764e56SPeter Brune 
826f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx);
827f190f2fcSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx);
8289e764e56SPeter Brune 
8299bd66eb0SPeter Brune /* set the functions for VI-specific line search operations */
8309bd66eb0SPeter Brune 
831*d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn *, SNESLineSearchVINormFn *, SNESLineSearchVIDirDerivFn *);
832*d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn **, SNESLineSearchVINormFn **, SNESLineSearchVIDirDerivFn **);
8339bd66eb0SPeter Brune 
8349e764e56SPeter Brune /* pointers to the associated SNES in order to be able to get the function evaluation out */
835014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch, SNES);
836014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch, SNES *);
8379e764e56SPeter Brune 
8389e764e56SPeter Brune /* set and get the parameters and vectors */
839014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
840014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscInt);
8419e764e56SPeter Brune 
842014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch, Vec, Vec, PetscBool *, void *);
84386d74e61SPeter Brune 
844014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch, PetscReal *);
845014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch, PetscReal);
8469e764e56SPeter Brune 
847014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch, PetscReal *);
848014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch, PetscReal);
8499e764e56SPeter Brune 
850014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch, PetscInt *order);
851014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch, PetscInt order);
85259405d5eSPeter Brune 
853422a814eSBarry Smith /*E
854dc4c0fb0SBarry Smith     SNESLineSearchReason - indication if the line search has succeeded or failed and why
855422a814eSBarry Smith 
856dc4c0fb0SBarry Smith   Values:
857dc4c0fb0SBarry Smith +  `SNES_LINESEARCH_SUCCEEDED`       - the line search succeeded
858dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_NANORINF` - a not a number of infinity appeared in the computions
859dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_DOMAIN`   - the function was evaluated outside of its domain, see `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
860dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_REDUCT`   - the linear search failed to get the requested decrease in its norm or objective
861dc4c0fb0SBarry Smith .  `SNES_LINESEARCH_FAILED_USER`     - used by `SNESLINESEARCHNLEQERR` to indicate the user changed the search direction inappropriately
862dc4c0fb0SBarry Smith -  `SNES_LINESEARCH_FAILED_FUNCTION` - indicates the maximum number of function evaluations allowed has been surpassed, `SNESConvergedReason` is also
863dc4c0fb0SBarry Smith                                        set to `SNES_DIVERGED_FUNCTION_COUNT`
864422a814eSBarry Smith 
86516a05f60SBarry Smith    Level: intermediate
86616a05f60SBarry Smith 
867dc4c0fb0SBarry Smith    Developer Note:
868dc4c0fb0SBarry Smith    Some of these reasons overlap with values of `SNESConvergedReason`
869422a814eSBarry Smith 
8701cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`,
871dc4c0fb0SBarry Smith           `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
872422a814eSBarry Smith E*/
8739371c9d4SSatish Balay typedef enum {
8749371c9d4SSatish Balay   SNES_LINESEARCH_SUCCEEDED,
875422a814eSBarry Smith   SNES_LINESEARCH_FAILED_NANORINF,
876e9b602ebSSatish Balay   SNES_LINESEARCH_FAILED_DOMAIN,
877d5b43468SJose E. Roman   SNES_LINESEARCH_FAILED_REDUCT, /* INSUFFICIENT REDUCTION */
878e9b602ebSSatish Balay   SNES_LINESEARCH_FAILED_USER,
8799371c9d4SSatish Balay   SNES_LINESEARCH_FAILED_FUNCTION
8809371c9d4SSatish Balay } SNESLineSearchReason;
881422a814eSBarry Smith 
882422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetReason(SNESLineSearch, SNESLineSearchReason *);
883422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetReason(SNESLineSearch, SNESLineSearchReason);
8849e764e56SPeter Brune 
885014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch, Vec *, Vec *, Vec *, Vec *, Vec *);
886014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch, Vec, Vec, Vec, Vec, Vec);
8879e764e56SPeter Brune 
888014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *);
889014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch, PetscReal, PetscReal, PetscReal);
890014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch);
891014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch, PetscBool);
8929e764e56SPeter Brune 
893dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitor(SNESLineSearch);
89449abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, void *), void *, PetscCtxDestroyFn *);
895d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch, const char[], const char[], const char[], PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *));
896dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch);
897dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch, PetscViewer);
898dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch, PetscViewer *);
899d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch, PetscViewerAndFormat *);
9009e764e56SPeter Brune 
901014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch, const char prefix[]);
902014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch, const char *prefix[]);
9039e764e56SPeter Brune 
9049e764e56SPeter Brune /* Shell interface functions */
9058434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellSetApply(SNESLineSearch, SNESLineSearchShellApplyFn, void *);
9068434afd1SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellGetApply(SNESLineSearch, SNESLineSearchShellApplyFn **, void **);
9079bcc50f1SBarry Smith 
9089bcc50f1SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellSetApply()", ) static inline PetscErrorCode SNESLineSearchShellSetUserFunc(SNESLineSearch ls, SNESLineSearchUserFunc f, void *ctx)
9099bcc50f1SBarry Smith {
9109bcc50f1SBarry Smith   return SNESLineSearchShellSetApply(ls, f, ctx);
9119bcc50f1SBarry Smith }
9129bcc50f1SBarry Smith 
9139bcc50f1SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellGetApply()", ) static inline PetscErrorCode SNESLineSearchShellGetUserFunc(SNESLineSearch ls, SNESLineSearchUserFunc *f, void **ctx)
9149bcc50f1SBarry Smith {
9159bcc50f1SBarry Smith   return SNESLineSearchShellGetApply(ls, f, ctx);
9169bcc50f1SBarry Smith }
9179e764e56SPeter Brune 
9182f4102e2SPeter Brune /* BT interface functions */
919014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch, PetscReal);
920014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch, PetscReal *);
9212f4102e2SPeter Brune 
9229e764e56SPeter Brune /*register line search types */
923bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchRegister(const char[], PetscErrorCode (*)(SNESLineSearch));
9249e764e56SPeter Brune 
925720c9a41SShri Abhyankar /* Routines for VI solver */
926014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetVariableBounds(SNES, Vec, Vec);
927cf836535SBlaise Bourdin PETSC_EXTERN PetscErrorCode SNESVIGetVariableBounds(SNES, Vec *, Vec *);
928014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
929014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetInactiveSet(SNES, IS *);
930014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetActiveSetIS(SNES, Vec, Vec, IS *);
931014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES, Vec, Vec, PetscReal *);
932*d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFtY(SNES, Vec, Vec, Vec, PetscScalar *);
933014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetRedundancyCheck(SNES, PetscErrorCode (*)(SNES, IS, IS *, void *), void *);
934f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeMeritFunction(Vec, PetscReal *, PetscReal *);
935f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeFunction(SNES, Vec, Vec, void *);
936f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMSetVI(DM, IS);
937f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMDestroyVI(DM);
938f5ea5bd2SShri Abhyankar 
939014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTestLocalMin(SNES);
940da9b6338SBarry Smith 
941eef9c623SLois Curfman McInnes /* Should this routine be private? */
942d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian(SNES, Vec, Mat, Mat);
943e885f1abSBarry Smith PETSC_EXTERN PetscErrorCode SNESTestJacobian(SNES);
944494a190aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESTestFunction(SNES);
945ddbbdb52SLois Curfman McInnes 
946014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetDM(SNES, DM);
947014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetDM(SNES, DM *);
948be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPC(SNES, SNES);
949be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPC(SNES, SNES *);
9503ad1a0b9SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESHasNPC(SNES, PetscBool *);
951be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESApplyNPC(SNES, Vec, Vec, Vec);
952be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCFunction(SNES, Vec, PetscReal *);
953be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeFunctionDefaultNPC(SNES, Vec, Vec);
954be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPCSide(SNES, PCSide);
955be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCSide(SNES, PCSide *);
9567601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLineSearch(SNES, SNESLineSearch);
9577601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLineSearch(SNES, SNESLineSearch *);
9586c699258SBarry Smith 
959edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESGetLineSearch()", ) static inline PetscErrorCode SNESGetSNESLineSearch(SNES snes, SNESLineSearch *ls)
960d71ae5a4SJacob Faibussowitsch {
9619371c9d4SSatish Balay   return SNESGetLineSearch(snes, ls);
9629371c9d4SSatish Balay }
963edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESSetLineSearch()", ) static inline PetscErrorCode SNESSetSNESLineSearch(SNES snes, SNESLineSearch ls)
964d71ae5a4SJacob Faibussowitsch {
9659371c9d4SSatish Balay   return SNESSetLineSearch(snes, ls);
9669371c9d4SSatish Balay }
9678b7b3213SJed Brown 
968014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetUpMatrices(SNES);
9698434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunction(DM, SNESFunctionFn *, void *);
9708434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetFunction(DM, SNESFunctionFn **, void **);
97149abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionContextDestroy(DM, PetscCtxDestroyFn *);
9728434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetMFFunction(DM, SNESFunctionFn *, void *);
9738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetNGS(DM, SNESNGSFn *, void *);
9748434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetNGS(DM, SNESNGSFn **, void **);
9758434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobian(DM, SNESJacobianFn *, void *);
9768434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetJacobian(DM, SNESJacobianFn **, void **);
97749abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianContextDestroy(DM, PetscCtxDestroyFn *);
9788434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetPicard(DM, SNESFunctionFn *, SNESJacobianFn *, void *);
9798434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetPicard(DM, SNESFunctionFn **, SNESJacobianFn **, void **);
9808434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetObjective(DM, SNESObjectiveFn *, void *);
9818434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetObjective(DM, SNESObjectiveFn **, void **);
9824dd50a75SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMSNES(DM, DM);
9836cab3a1bSJed Brown 
9848434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDASNESFunctionFn)(DMDALocalInfo *, void *, void *, void *);
9858434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDASNESJacobianFn)(DMDALocalInfo *, void *, Mat, Mat, void *);
9868434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDASNESObjectiveFn)(DMDALocalInfo *, void *, PetscReal *, void *);
987c9d099b5SPeter Brune 
9888434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDASNESFunctionVecFn)(DMDALocalInfo *, Vec, Vec, void *);
9898434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDASNESJacobianVecFn)(DMDALocalInfo *, Vec, Mat, Mat, void *);
9908434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDASNESObjectiveVecFn)(DMDALocalInfo *, Vec, PetscReal *, void *);
9915505f7afSJunchao Zhang 
9928434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocal(DM, InsertMode, DMDASNESFunctionFn *, void *);
9938434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocal(DM, DMDASNESJacobianFn *, void *);
9948434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocal(DM, DMDASNESObjectiveFn *, void *);
9958434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetPicardLocal(DM, InsertMode, DMDASNESFunctionFn *, DMDASNESJacobianFn, void *);
9966cab3a1bSJed Brown 
9978434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocalVec(DM, InsertMode, DMDASNESFunctionVecFn *, void *);
9988434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocalVec(DM, DMDASNESJacobianVecFn *, void *);
9998434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocalVec(DM, DMDASNESObjectiveVecFn *, void *);
10005505f7afSJunchao Zhang 
1001bdd6f66aSToby Isaac PETSC_EXTERN PetscErrorCode DMSNESSetBoundaryLocal(DM, PetscErrorCode (*)(DM, Vec, void *), void *);
10026493148fSStefano Zampini PETSC_EXTERN PetscErrorCode DMSNESSetObjectiveLocal(DM, PetscErrorCode (*)(DM, Vec, PetscReal *, void *), void *);
1003ff35dfedSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionLocal(DM, PetscErrorCode (*)(DM, Vec, Vec, void *), void *);
1004d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianLocal(DM, PetscErrorCode (*)(DM, Vec, Mat, Mat, void *), void *);
100528d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetBoundaryLocal(DM, PetscErrorCode (**)(DM, Vec, void *), void **);
10066493148fSStefano Zampini PETSC_EXTERN PetscErrorCode DMSNESGetObjectiveLocal(DM, PetscErrorCode (**)(DM, Vec, PetscReal *, void *), void **);
100728d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetFunctionLocal(DM, PetscErrorCode (**)(DM, Vec, Vec, void *), void **);
100828d58a37SPierre Jolivet PETSC_EXTERN PetscErrorCode DMSNESGetJacobianLocal(DM, PetscErrorCode (**)(DM, Vec, Mat, Mat, void *), void **);
1009ff35dfedSBarry Smith 
1010b665c654SMatthew G Knepley /* Routines for Multiblock solver */
1011014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetFields(SNES, const char[], PetscInt, const PetscInt *);
1012014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetIS(SNES, const char[], IS);
1013014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
1014014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);
10154bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMultiblockGetSubSNES(SNES, PetscInt *, SNES *[]);
1016aeed3662SMatthew G Knepley 
101737e1895aSJed Brown /*J
101887497f52SBarry Smith    SNESMSType - String with the name of a PETSc `SNESMS` method.
101937e1895aSJed Brown 
102037e1895aSJed Brown    Level: intermediate
102137e1895aSJed Brown 
10221cc06b55SBarry Smith .seealso: [](ch_snes), `SNESMS`, `SNESMSGetType()`, `SNESMSSetType()`, `SNES`
102337e1895aSJed Brown J*/
102419fd82e9SBarry Smith typedef const char *SNESMSType;
102537e1895aSJed Brown #define SNESMSM62       "m62"
102637e1895aSJed Brown #define SNESMSEULER     "euler"
1027a97cb6bcSJed Brown #define SNESMSJAMESON83 "jameson83"
10283847c725SLisandro Dalcin #define SNESMSVLTP11    "vltp11"
1029a97cb6bcSJed Brown #define SNESMSVLTP21    "vltp21"
1030a97cb6bcSJed Brown #define SNESMSVLTP31    "vltp31"
1031a97cb6bcSJed Brown #define SNESMSVLTP41    "vltp41"
1032a97cb6bcSJed Brown #define SNESMSVLTP51    "vltp51"
1033a97cb6bcSJed Brown #define SNESMSVLTP61    "vltp61"
103437e1895aSJed Brown 
103519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSRegister(SNESMSType, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[]);
10364bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMSRegisterAll(void);
103757715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetType(SNES, SNESMSType *);
103819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSSetType(SNES, SNESMSType);
103957715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetDamping(SNES, PetscReal *);
104057715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSSetDamping(SNES, PetscReal);
1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSFinalizePackage(void);
1042607a6623SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSInitializePackage(void);
1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSRegisterDestroy(void);
104437e1895aSJed Brown 
1045ceaaa498SBarry Smith /*MC
1046ceaaa498SBarry Smith    SNESNGMRESRestartType - the restart approach used by `SNESNGMRES`
104713a62661SPeter Brune 
1048ceaaa498SBarry Smith   Values:
1049ceaaa498SBarry Smith +   `SNES_NGMRES_RESTART_NONE`       - never restart
1050ceaaa498SBarry Smith .   `SNES_NGMRES_RESTART_DIFFERENCE` - restart based upon difference criteria
1051ceaaa498SBarry Smith -   `SNES_NGMRES_RESTART_PERIODIC`   - restart after a fixed number of iterations
1052ceaaa498SBarry Smith 
1053ceaaa498SBarry Smith   Options Database Keys:
1054ceaaa498SBarry Smith + -snes_ngmres_restart_type <difference,periodic,none> - set the restart type
1055af27ebaaSBarry Smith - -snes_ngmres_restart <30>                            - sets the number of iterations before restart for periodic
1056ceaaa498SBarry Smith 
105795bd0b28SBarry Smith    Level: intermediate
105895bd0b28SBarry Smith 
1059ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1060ceaaa498SBarry Smith           `SNESNGMRESGetRestartType()`, `SNESNGMRESSelectType`
1061ceaaa498SBarry Smith M*/
106213a62661SPeter Brune typedef enum {
106313a62661SPeter Brune   SNES_NGMRES_RESTART_NONE       = 0,
106413a62661SPeter Brune   SNES_NGMRES_RESTART_PERIODIC   = 1,
10659371c9d4SSatish Balay   SNES_NGMRES_RESTART_DIFFERENCE = 2
10669371c9d4SSatish Balay } SNESNGMRESRestartType;
10676a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];
106813a62661SPeter Brune 
1069ceaaa498SBarry Smith /*MC
1070ceaaa498SBarry Smith    SNESNGMRESSelectType - the approach used by `SNESNGMRES` to determine how the candidate solution and
1071ceaaa498SBarry Smith   combined solution are used to create the next iterate.
1072ceaaa498SBarry Smith 
1073ceaaa498SBarry Smith    Values:
1074ceaaa498SBarry Smith +   `SNES_NGMRES_SELECT_NONE`       - choose the combined solution all the time
1075ceaaa498SBarry Smith .   `SNES_NGMRES_SELECT_DIFFERENCE` - choose based upon the selection criteria
1076ceaaa498SBarry Smith -   `SNES_NGMRES_SELECT_LINESEARCH` - choose based upon line search combination
1077ceaaa498SBarry Smith 
1078ceaaa498SBarry Smith   Options Database Key:
1079ceaaa498SBarry Smith . -snes_ngmres_select_type<difference,none,linesearch> - select type
1080ceaaa498SBarry Smith 
108195bd0b28SBarry Smith    Level: intermediate
108295bd0b28SBarry Smith 
1083ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1084ceaaa498SBarry Smith           `SNESNGMRESGetRestartType()`, `SNESNGMRESRestartType`
1085ceaaa498SBarry Smith M*/
108613a62661SPeter Brune typedef enum {
108713a62661SPeter Brune   SNES_NGMRES_SELECT_NONE       = 0,
108813a62661SPeter Brune   SNES_NGMRES_SELECT_DIFFERENCE = 1,
10899371c9d4SSatish Balay   SNES_NGMRES_SELECT_LINESEARCH = 2
10909371c9d4SSatish Balay } SNESNGMRESSelectType;
10916a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESSelectTypes[];
109213a62661SPeter Brune 
1093014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartType(SNES, SNESNGMRESRestartType);
1094014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetSelectType(SNES, SNESNGMRESSelectType);
109523b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartFmRise(SNES, PetscBool);
109623b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESGetRestartFmRise(SNES, PetscBool *);
109713a62661SPeter Brune 
1098ceaaa498SBarry Smith /*MC
1099ceaaa498SBarry Smith    SNESNCGType - the conjugate update approach for `SNESNCG`
11000a844d1aSPeter Brune 
1101ceaaa498SBarry Smith    Values:
1102ceaaa498SBarry Smith +   `SNES_NCG_FR`  - Fletcher-Reeves update
1103ceaaa498SBarry Smith .   `SNES_NCG_PRP` - Polak-Ribiere-Polyak update, the default and the only one that tolerates generalized search directions
1104ceaaa498SBarry Smith .   `SNES_NCG_HS`  - Hestenes-Steifel update
1105ceaaa498SBarry Smith .   `SNES_NCG_DY`  - Dai-Yuan update
1106ceaaa498SBarry Smith -   `SNES_NCG_CD`  - Conjugate Descent update
1107ceaaa498SBarry Smith 
1108ceaaa498SBarry Smith   Options Database Key:
1109ceaaa498SBarry Smith . -snes_ncg_type<fr,prp,hs,dy,cd> - select type
1110ceaaa498SBarry Smith 
111195bd0b28SBarry Smith    Level: intermediate
111295bd0b28SBarry Smith 
1113ceaaa498SBarry Smith .seealso: `SNES, `SNESNCG`, `SNESNCGSetType()`
1114ceaaa498SBarry Smith M*/
11150a844d1aSPeter Brune typedef enum {
11160de8b71eSPeter Brune   SNES_NCG_FR  = 0,
11170de8b71eSPeter Brune   SNES_NCG_PRP = 1,
11180de8b71eSPeter Brune   SNES_NCG_HS  = 2,
11190de8b71eSPeter Brune   SNES_NCG_DY  = 3,
11209371c9d4SSatish Balay   SNES_NCG_CD  = 4
11219371c9d4SSatish Balay } SNESNCGType;
11226a6fc655SJed Brown PETSC_EXTERN const char *const SNESNCGTypes[];
11230a844d1aSPeter Brune 
1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNCGSetType(SNES, SNESNCGType);
112513a62661SPeter Brune 
1126ceaaa498SBarry Smith /*MC
1127ceaaa498SBarry Smith    SNESQNScaleType - the scaling type used by `SNESQN`
1128ceaaa498SBarry Smith 
1129ceaaa498SBarry Smith    Values:
1130ceaaa498SBarry Smith +   `SNES_QN_SCALE_NONE`     - don't scale the problem
1131ceaaa498SBarry Smith .   `SNES_QN_SCALE_SCALAR`   - use Shanno scaling
1132ceaaa498SBarry Smith .   `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
1133ceaaa498SBarry Smith -   `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with `SNESSetJacobian()`
1134ceaaa498SBarry Smith                                computed at the first iteration of `SNESQN` and at ever restart.
1135ceaaa498SBarry Smith 
1136ceaaa498SBarry Smith     Options Database Key:
1137ceaaa498SBarry Smith . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
1138ceaaa498SBarry Smith 
113995bd0b28SBarry Smith    Level: intermediate
114095bd0b28SBarry Smith 
1141ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNRestartType`
1142ceaaa498SBarry Smith M*/
11439371c9d4SSatish Balay typedef enum {
11449371c9d4SSatish Balay   SNES_QN_SCALE_DEFAULT  = 0,
11451efc8c45SPeter Brune   SNES_QN_SCALE_NONE     = 1,
114692f76d53SAlp Dener   SNES_QN_SCALE_SCALAR   = 2,
114792f76d53SAlp Dener   SNES_QN_SCALE_DIAGONAL = 3,
11489371c9d4SSatish Balay   SNES_QN_SCALE_JACOBIAN = 4
11499371c9d4SSatish Balay } SNESQNScaleType;
11506a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNScaleTypes[];
1151ceaaa498SBarry Smith 
1152ceaaa498SBarry Smith /*MC
1153ceaaa498SBarry Smith    SNESQNRestartType - the restart approached used by `SNESQN`
1154ceaaa498SBarry Smith 
1155ceaaa498SBarry Smith    Values:
1156ceaaa498SBarry Smith +   `SNES_QN_RESTART_NONE`     - never restart
1157ceaaa498SBarry Smith .   `SNES_QN_RESTART_POWELL`   - restart based upon descent criteria
1158ceaaa498SBarry Smith -   `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations
1159ceaaa498SBarry Smith 
1160ceaaa498SBarry Smith   Options Database Keys:
1161ceaaa498SBarry Smith + -snes_qn_restart_type <powell,periodic,none> - set the restart type
1162ceaaa498SBarry Smith - -snes_qn_m <m>                               - sets the number of stored updates and the restart period for periodic
1163ceaaa498SBarry Smith 
116495bd0b28SBarry Smith    Level: intermediate
116595bd0b28SBarry Smith 
1166ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNScaleType`
1167ceaaa498SBarry Smith M*/
11689371c9d4SSatish Balay typedef enum {
11699371c9d4SSatish Balay   SNES_QN_RESTART_DEFAULT  = 0,
11701efc8c45SPeter Brune   SNES_QN_RESTART_NONE     = 1,
11711efc8c45SPeter Brune   SNES_QN_RESTART_POWELL   = 2,
11729371c9d4SSatish Balay   SNES_QN_RESTART_PERIODIC = 3
11739371c9d4SSatish Balay } SNESQNRestartType;
11746a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNRestartTypes[];
1175ceaaa498SBarry Smith 
1176ceaaa498SBarry Smith /*MC
1177ceaaa498SBarry Smith    SNESQNType - the type used by `SNESQN`
1178ceaaa498SBarry Smith 
1179ceaaa498SBarry Smith   Values:
1180ceaaa498SBarry Smith +   `SNES_QN_LBFGS`      - LBFGS variant
1181ceaaa498SBarry Smith .   `SNES_QN_BROYDEN`    - Broyden variant
1182ceaaa498SBarry Smith -   `SNES_QN_BADBROYDEN` - Bad Broyden variant
1183ceaaa498SBarry Smith 
1184ceaaa498SBarry Smith   Options Database Key:
1185ceaaa498SBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
1186ceaaa498SBarry Smith 
118795bd0b28SBarry Smith    Level: intermediate
118895bd0b28SBarry Smith 
1189ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNSetType()`, `SNESQNScaleType`, `SNESQNRestartType`, `SNESQNSetRestartType()`
1190ceaaa498SBarry Smith M*/
11919371c9d4SSatish Balay typedef enum {
11929371c9d4SSatish Balay   SNES_QN_LBFGS      = 0,
11931efc8c45SPeter Brune   SNES_QN_BROYDEN    = 1,
11941efc8c45SPeter Brune   SNES_QN_BADBROYDEN = 2
11951efc8c45SPeter Brune } SNESQNType;
11961efc8c45SPeter Brune PETSC_EXTERN const char *const SNESQNTypes[];
11970c777b0cSPeter Brune 
11981efc8c45SPeter Brune PETSC_EXTERN PetscErrorCode SNESQNSetType(SNES, SNESQNType);
1199014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetScaleType(SNES, SNESQNScaleType);
1200014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetRestartType(SNES, SNESQNRestartType);
12010c777b0cSPeter Brune 
1202e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetType(SNES, PCASMType *);
1203e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetType(SNES, PCASMType);
120476857b2aSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomains(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **);
120576857b2aSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetSubdomains(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *);
1206610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetDamping(SNES, PetscReal);
1207610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetDamping(SNES, PetscReal *);
120876857b2aSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomainVecs(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **);
1209d728fb7dSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES, PetscBool);
1210d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetSNES(SNES, PetscInt, SNES *);
1211d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetNumber(SNES, PetscInt *);
1212f10b3e88SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMSetWeight(SNES, Vec);
12130c777b0cSPeter Brune 
1214ad05f2b7SBarry Smith /*E
1215ad05f2b7SBarry Smith   SNESCompositeType - Determines how two or more preconditioners are composed with the `SNESType` of `SNESCOMPOSITE`
1216ad05f2b7SBarry Smith 
1217ad05f2b7SBarry Smith   Values:
1218ad05f2b7SBarry Smith + `SNES_COMPOSITE_ADDITIVE`        - results from application of all preconditioners are added together
1219ad05f2b7SBarry Smith . `SNES_COMPOSITE_MULTIPLICATIVE`  - preconditioners are applied sequentially to the residual freshly
1220ad05f2b7SBarry Smith                                      computed after the previous preconditioner application
1221ad05f2b7SBarry Smith - `SNES_COMPOSITE_ADDITIVEOPTIMAL` - uses a linear combination of the solutions obtained with each preconditioner that approximately minimize the function
1222ad05f2b7SBarry Smith                                      value at the new iteration.
1223ad05f2b7SBarry Smith 
1224ad05f2b7SBarry Smith    Level: beginner
1225ad05f2b7SBarry Smith 
1226ad05f2b7SBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `PCCompositeType`
1227ad05f2b7SBarry Smith E*/
12289371c9d4SSatish Balay typedef enum {
12299371c9d4SSatish Balay   SNES_COMPOSITE_ADDITIVE,
12309371c9d4SSatish Balay   SNES_COMPOSITE_MULTIPLICATIVE,
12319371c9d4SSatish Balay   SNES_COMPOSITE_ADDITIVEOPTIMAL
12329371c9d4SSatish Balay } SNESCompositeType;
1233eed5f15bSPeter Brune PETSC_EXTERN const char *const SNESCompositeTypes[];
1234eed5f15bSPeter Brune 
1235eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetType(SNES, SNESCompositeType);
1236eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeAddSNES(SNES, SNESType);
1237eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeGetSNES(SNES, PetscInt, SNES *);
1238a6b47ab3SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESCompositeGetNumber(SNES, PetscInt *);
12398f626970SPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetDamping(SNES, PetscInt, PetscReal);
1240eed5f15bSPeter Brune 
124168e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetDiscretisationInfo(SNES, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
12424d04e9f1SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetComputeOperator(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
124339fd2e8aSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetComputeFunction(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, void *), void *);
124468e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetConstructType(SNES, PCPatchConstructType, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, void *), void *);
124568e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetCellNumbering(SNES, PetscSection);
1246561742edSMatthew G. Knepley 
1247a055c8caSBarry Smith /*E
1248a055c8caSBarry Smith     SNESFASType - Determines the type of nonlinear multigrid method that is run.
1249a055c8caSBarry Smith 
1250a055c8caSBarry Smith    Values:
125187497f52SBarry Smith +  `SNES_FAS_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `SNESFASSetCycles()`
125287497f52SBarry Smith .  `SNES_FAS_ADDITIVE`                 - additive FAS cycle
125387497f52SBarry Smith .  `SNES_FAS_FULL`                     - full FAS cycle
125487497f52SBarry Smith -  `SNES_FAS_KASKADE`                  - Kaskade FAS cycle
1255a055c8caSBarry Smith 
125616a05f60SBarry Smith    Level: beginner
125716a05f60SBarry Smith 
12581cc06b55SBarry Smith .seealso: [](ch_snes), `SNESFAS`, `PCMGSetType()`, `PCMGType`
1259a055c8caSBarry Smith E*/
12609371c9d4SSatish Balay typedef enum {
12619371c9d4SSatish Balay   SNES_FAS_MULTIPLICATIVE,
12629371c9d4SSatish Balay   SNES_FAS_ADDITIVE,
12639371c9d4SSatish Balay   SNES_FAS_FULL,
12649371c9d4SSatish Balay   SNES_FAS_KASKADE
12659371c9d4SSatish Balay } SNESFASType;
1266a055c8caSBarry Smith PETSC_EXTERN const char *const SNESFASTypes[];
1267a055c8caSBarry Smith 
1268a055c8caSBarry Smith /* called on the finest level FAS instance*/
1269a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetType(SNES, SNESFASType);
1270a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetType(SNES, SNESFASType *);
1271a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLevels(SNES, PetscInt, MPI_Comm *);
1272a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetLevels(SNES, PetscInt *);
1273a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCycleSNES(SNES, PetscInt, SNES *);
1274a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothUp(SNES, PetscInt);
1275a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothDown(SNES, PetscInt);
1276a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetCycles(SNES, PetscInt);
1277d142ab34SLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESFASSetMonitor(SNES, PetscViewerAndFormat *, PetscBool);
1278a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLog(SNES, PetscBool);
1279a055c8caSBarry Smith 
1280a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetGalerkin(SNES, PetscBool);
1281a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetGalerkin(SNES, PetscBool *);
128225acbd8eSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESFASGalerkinFunctionDefault(SNES, Vec, Vec, void *);
1283a055c8caSBarry Smith 
1284a055c8caSBarry Smith /* called on any level -- "Cycle" FAS instance */
1285a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmoother(SNES, SNES *);
1286a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherUp(SNES, SNES *);
1287a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherDown(SNES, SNES *);
1288a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetCorrection(SNES, SNES *);
1289a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInterpolation(SNES, Mat *);
1290a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRestriction(SNES, Mat *);
1291a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInjection(SNES, Mat *);
1292a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRScale(SNES, Vec *);
1293a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleSetCycles(SNES, PetscInt);
1294a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleIsFine(SNES, PetscBool *);
1295a055c8caSBarry Smith 
1296a055c8caSBarry Smith /* called on the (outer) finest level FAS to set/get parameters on any level instance */
1297a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInterpolation(SNES, PetscInt, Mat);
1298a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInterpolation(SNES, PetscInt, Mat *);
1299a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRestriction(SNES, PetscInt, Mat);
1300a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRestriction(SNES, PetscInt, Mat *);
1301a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInjection(SNES, PetscInt, Mat);
1302a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInjection(SNES, PetscInt, Mat *);
1303a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRScale(SNES, PetscInt, Vec);
1304a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRScale(SNES, PetscInt, Vec *);
1305020631bcSPeter Brune PETSC_EXTERN PetscErrorCode SNESFASSetContinuation(SNES, PetscBool);
1306a055c8caSBarry Smith 
1307a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmoother(SNES, PetscInt, SNES *);
1308a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherUp(SNES, PetscInt, SNES *);
1309a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherDown(SNES, PetscInt, SNES *);
1310a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCoarseSolve(SNES, SNES *);
1311a055c8caSBarry Smith 
1312a055c8caSBarry Smith /* parameters for full FAS */
1313a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASFullSetDownSweep(SNES, PetscBool);
1314a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCreateCoarseVec(SNES, Vec *);
1315a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASRestrict(SNES, Vec, Vec);
13166dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullSetTotal(SNES, PetscBool);
13176dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullGetTotal(SNES, PetscBool *);
1318a055c8caSBarry Smith 
13192eace19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetSNESVariableBounds(DM, SNES);
1320f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckDiscretization(SNES, DM, PetscReal, Vec, PetscReal, PetscReal[]);
1321f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckResidual(SNES, DM, Vec, PetscReal, PetscReal *);
1322f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckJacobian(SNES, DM, Vec, PetscReal, PetscBool *, PetscReal *);
13237f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckFromOptions(SNES, Vec);
13248e3a2eefSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESComputeJacobianAction(DM, Vec, Vec, Vec, void *);
13258e3a2eefSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCreateJacobianMF(DM, Vec, void *, Mat *);
132697276fddSZach Atkins 
132797276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALSetFunction(SNES, SNESFunctionFn *, void *ctx);
132897276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALGetFunction(SNES, SNESFunctionFn **, void **ctx);
132997276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALComputeFunction(SNES, Vec, Vec);
133097276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALGetLoadParameter(SNES, PetscReal *);
133197276fddSZach Atkins 
133297276fddSZach Atkins /*MC
133397276fddSZach Atkins    SNESNewtonALCorrectionType - the approach used by `SNESNEWTONAL` to determine
133497276fddSZach Atkins    the correction to the current increment. While the exact correction satisfies
133597276fddSZach Atkins    the constraint surface at every iteration, it also requires solving a quadratic
133697276fddSZach Atkins    equation which may not have real roots. Conversely, the normal correction is more
133797276fddSZach Atkins    efficient and always yields a real correction and is the default.
133897276fddSZach Atkins 
133997276fddSZach Atkins    Values:
134097276fddSZach Atkins +   `SNES_NEWTONAL_CORRECTION_EXACT` - choose the correction which exactly satisfies the constraint
1341d7c1f440SPierre Jolivet -   `SNES_NEWTONAL_CORRECTION_NORMAL` - choose the correction in the updated normal hyper-surface to the constraint surface
134297276fddSZach Atkins 
134397276fddSZach Atkins    Options Database Key:
134497276fddSZach Atkins . -snes_newtonal_correction_type <exact> - select type from <exact,normal>
134597276fddSZach Atkins 
134697276fddSZach Atkins    Level: intermediate
134797276fddSZach Atkins 
134897276fddSZach Atkins .seealso: `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetCorrectionType()`
134997276fddSZach Atkins M*/
135097276fddSZach Atkins typedef enum {
135197276fddSZach Atkins   SNES_NEWTONAL_CORRECTION_EXACT  = 0,
135297276fddSZach Atkins   SNES_NEWTONAL_CORRECTION_NORMAL = 1,
135397276fddSZach Atkins } SNESNewtonALCorrectionType;
135497276fddSZach Atkins PETSC_EXTERN const char *const SNESNewtonALCorrectionTypes[];
135597276fddSZach Atkins 
135697276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALSetCorrectionType(SNES, SNESNewtonALCorrectionType);
1357