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);
64*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorSAWs(SNES, PetscInt, PetscReal, PetscCtx);
658a70d858SHong Zhang PETSC_EXTERN PetscErrorCode SNESMonitorSAWsCreate(SNES, void **);
66*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESMonitorSAWsDestroy(PetscCtxRt);
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
SNESReasonView(SNES snes,PetscViewer v)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 }
SNESReasonViewFromOptions(SNES snes)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
169*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRSetPreCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, PetscBool *, void *), PetscCtx);
170*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRGetPreCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, PetscBool *, void *), PetscCtxRt);
171*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRSetPostCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), PetscCtx);
172*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRGetPostCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), PetscCtxRt);
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 *);
207*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, PetscBool *, PetscCtx), PetscCtx);
208*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, PetscBool *, PetscCtx), PetscCtxRt);
209*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES, PetscErrorCode (*)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, PetscCtx), PetscCtx);
210*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES, PetscErrorCode (**)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, PetscCtx), PetscCtxRt);
21141ba4c6cSHeeho Park
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetNonlinearStepFailures(SNES, PetscInt *);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetMaxNonlinearStepFailures(SNES, PetscInt);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetMaxNonlinearStepFailures(SNES, PetscInt *);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetNumberFunctionEvals(SNES, PetscInt *);
216b850b91aSLisandro Dalcin
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLagPreconditioner(SNES, PetscInt);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLagPreconditioner(SNES, PetscInt *);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLagJacobian(SNES, PetscInt);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLagJacobian(SNES, PetscInt *);
22137ec4e1aSPeter Brune PETSC_EXTERN PetscErrorCode SNESSetLagPreconditionerPersists(SNES, PetscBool);
22237ec4e1aSPeter Brune PETSC_EXTERN PetscErrorCode SNESSetLagJacobianPersists(SNES, PetscBool);
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetGridSequence(SNES, PetscInt);
224fa19ca70SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetGridSequence(SNES, PetscInt *);
225a8054027SBarry Smith
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLinearSolveIterations(SNES, PetscInt *);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLinearSolveFailures(SNES, PetscInt *);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetMaxLinearSolveFailures(SNES, PetscInt);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetMaxLinearSolveFailures(SNES, PetscInt *);
230971e163fSPeter Brune PETSC_EXTERN PetscErrorCode SNESSetCountersReset(SNES, PetscBool);
23112b1dd1aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESResetCounters(SNES);
2323d4c4710SBarry Smith
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPSetUseEW(SNES, PetscBool);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPGetUseEW(SNES, PetscBool *);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPSetParametersEW(SNES, PetscInt, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESKSPGetParametersEW(SNES, PetscInt *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
237eafb4bcbSBarry Smith
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMonitorLGRange(SNES, PetscInt, PetscReal, void *);
239eafb4bcbSBarry Smith
240*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetApplicationContext(SNES, PetscCtx);
241*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetApplicationContext(SNES, PetscCtxRt);
242*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetComputeApplicationContext(SNES, PetscErrorCode (*)(SNES, PetscCtxRt), PetscCtxDestroyFn *);
243184914b5SBarry Smith
244014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESPythonSetType(SNES, const char[]);
245ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode SNESPythonGetType(SNES, const char *[]);
2461d6018f0SLisandro Dalcin
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetFunctionDomainError(SNES);
24876c63389SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetObjectiveDomainError(SNES);
249cc6b0f04SFande Kong PETSC_EXTERN PetscErrorCode SNESSetJacobianDomainError(SNES);
250b351a90bSFande Kong PETSC_EXTERN PetscErrorCode SNESSetCheckJacobianDomainError(SNES, PetscBool);
2518383d7d7SFande Kong PETSC_EXTERN PetscErrorCode SNESGetCheckJacobianDomainError(SNES, PetscBool *);
2526a388c36SPeter Brune
253edd03b47SJacob Faibussowitsch #define SNES_CONVERGED_TR_DELTA_DEPRECATED SNES_CONVERGED_TR_DELTA PETSC_DEPRECATED_ENUM(3, 12, 0, "SNES_DIVERGED_TR_DELTA", )
2545664c22fSDavid Wells #define SNES_DIVERGED_FNORM_NAN_DEPRECATED SNES_DIVERGED_FNORM_NAN PETSC_DEPRECATED_ENUM(3, 25, 0, "SNES_DIVERGED_FUNCTION_NANORINF", )
255435da068SBarry Smith /*E
256dc4c0fb0SBarry Smith SNESConvergedReason - reason a `SNESSolve()` was determined to have converged or diverged
257435da068SBarry Smith
258dc4c0fb0SBarry Smith Values:
25976c63389SBarry Smith + `SNES_CONVERGED_FNORM_ABS` - $ ||F|| \le abstol $
26076c63389SBarry Smith . `SNES_CONVERGED_FNORM_RELATIVE` - $ ||F|| <= rtol*||F(x_0))|| $ where $x_0 $ is the initial guess
26176c63389SBarry Smith . `SNES_CONVERGED_SNORM_RELATIVE` - The 2-norm of the last step $ \le stol * ||x|| $ where $ x $ is the current solution
262379ef8d2SAlexander . `SNES_CONVERGED_USER` - The user has indicated convergence for an arbitrary reason
263dc4c0fb0SBarry Smith . `SNES_DIVERGED_FUNCTION_COUNT` - The user provided function has been called more times than the maximum set in `SNESSetTolerances()`
264dc4c0fb0SBarry Smith . `SNES_DIVERGED_DTOL` - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
26576c63389SBarry Smith . `SNES_DIVERGED_FUNCTION_NANORINF` - the 2-norm of the current function evaluation is not-a-number (NaN) or infinity, (this
26676c63389SBarry Smith is usually caused by a division of 0 by 0) and the solver could not recover from this (by, for example, cutting the step size)
26776c63389SBarry Smith . `SNES_DIVERGED_OBJECTIVE_NANORINF` - the object function evaluation is not-a-number (NaN) or infinity, (this
26876c63389SBarry Smith is usually caused by a division of 0 by 0) and the solver could not recover from this (by, for example, cutting the step size)
26976c63389SBarry Smith . `SNES_DIVERGED_FUNCTION_DOMAIN` - the function evaluation occurred outside the function's domain (function callback provided by
27076c63389SBarry Smith `SNESSetFunction()` called `SNESSetObjectiveDomainError()`) and the solver could not recover from this (by, for example, cutting the step size)
27176c63389SBarry Smith . `SNES_DIVERGED_OBJECTIVE_DOMAIN` - the object function evaluation occurred outside the function's domain (function callback provided by
27276c63389SBarry Smith `SNESSetObjective()` called `SNESSetObjectiveDomainError()`) and the solver could not recover from this (by, for example, cutting the step size)
27376c63389SBarry Smith . `SNES_DIVERGED_JACOBIAN_DOMAIN` - the Jacobian evaluation occurred outside the function's domain (function callback provided by
27476c63389SBarry Smith `SNESSetJacobian()` called `SNESSetJacobianDomainError()`)
275dc4c0fb0SBarry Smith . `SNES_DIVERGED_MAX_IT` - `SNESSolve()` has reached the maximum number of iterations requested
276dc4c0fb0SBarry Smith . `SNES_DIVERGED_LINE_SEARCH` - The line search has failed. This only occurs for `SNES` solvers that use a line search
277dc4c0fb0SBarry Smith . `SNES_DIVERGED_LOCAL_MIN` - the algorithm seems to have stagnated at a local minimum that is not zero.
278379ef8d2SAlexander - `SNES_CONVERGED_ITERATING - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
279f203c74bSBarry Smith
28016a05f60SBarry Smith Level: beginner
28116a05f60SBarry Smith
282dc4c0fb0SBarry Smith Notes:
2830b4b7b1cSBarry 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
2840b4b7b1cSBarry Smith (in this case we recommend
285dc4c0fb0SBarry Smith testing with `-pc_type lu` to eliminate the linear solver as the cause of the problem).
286dc4c0fb0SBarry Smith
2870b4b7b1cSBarry Smith `SNES_DIVERGED_LOCAL_MIN` can only occur when using a `SNES` solver that uses a line search (`SNESLineSearch`).
28876c63389SBarry Smith The line search wants to $ \min Q(\alpha) = 1/2 || F(x + \alpha s) ||^2_2 $ this occurs
28976c63389SBarry 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
29076c63389SBarry Smith $ Q'(\alpha) = -F(x)^T F'(x)^(-1)^T F'(x+\alpha s)F(x+\alpha s)$; when $\alpha = 0$
29176c63389SBarry Smith $Q'(0) = - ||F(x)||^2_2 $ which is always NEGATIVE if $F'(x)$ is invertible. This means the Newton
29276c63389SBarry Smith direction is a descent direction and the line search should succeed if $\alpha $ is small enough.
293f203c74bSBarry Smith
29476c63389SBarry Smith If $F'(x)$ is NOT invertible AND $F'(x)^T F(x) = 0 $ then $Q'(0) = 0 $ and the Newton direction
295f203c74bSBarry Smith is NOT a descent direction so the line search will fail. All one can do at this point
296f203c74bSBarry Smith is change the initial guess and try again.
297f203c74bSBarry Smith
298f203c74bSBarry Smith An alternative explanation: Newton's method can be regarded as replacing the function with
29976c63389SBarry Smith its linear approximation and minimizing the 2-norm of that. That is $F(x+s) \approx F(x) + F'(x)s$
30076c63389SBarry Smith so we minimize $ || F(x) + F'(x) s ||^2_2$ using Least Squares. If $F'(x)$ is invertible then
30176c63389SBarry 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
30276c63389SBarry Smith exists a nontrivial (that is $F'(x)s \ne 0$) solution to the equation and this direction is
30376c63389SBarry 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)
30476c63389SBarry Smith = - (F'(x)^T F(x)) [F'(x)^T F'(x)]^(-T) (F'(x)^T F(x))$. Since we are assuming $(F'(x)^T F(x)) \ne 0$
30576c63389SBarry Smith and $F'(x)^T F'(x)$ has no negative eigenvalues $Q'(0) < 0$ so $s$ is a descent direction and the line
30676c63389SBarry Smith search should succeed for small enough $\alpha$.
307f203c74bSBarry Smith
308f203c74bSBarry Smith Note that this RARELY happens in practice. Far more likely the linear system is not being solved
309f203c74bSBarry Smith (well enough?) or the Jacobian is wrong.
310f203c74bSBarry Smith
31187497f52SBarry Smith `SNES_DIVERGED_MAX_IT` means that the solver reached the maximum number of iterations without satisfying any
31287497f52SBarry Smith convergence criteria. `SNES_CONVERGED_ITS` means that `SNESConvergedSkip()` was chosen as the convergence test;
3134d0a8057SBarry Smith thus the usual convergence criteria have not been checked and may or may not be satisfied.
314f203c74bSBarry Smith
3151cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`, `SNESSetTolerances()`
316435da068SBarry Smith E*/
317184914b5SBarry Smith typedef enum { /* converged */
31801b82886SBarry Smith SNES_CONVERGED_FNORM_ABS = 2, /* ||F|| < atol */
31901b82886SBarry Smith SNES_CONVERGED_FNORM_RELATIVE = 3, /* ||F|| < rtol*||F_initial|| */
3205358d0d4SBarry Smith SNES_CONVERGED_SNORM_RELATIVE = 4, /* Newton computed step size small; || delta x || < stol || x || */
3213f149594SLisandro Dalcin SNES_CONVERGED_ITS = 5, /* maximum iterations reached */
322ce78bad3SBarry Smith SNES_BREAKOUT_INNER_ITER = 6, /* Flag to break out of inner loop after checking custom convergence, used in multi-phase flow when state changes */
323379ef8d2SAlexander SNES_CONVERGED_USER = 7, /* The user has indicated convergence for an arbitrary reason */
324184914b5SBarry Smith /* diverged */
32546a9e3ceSBarry Smith SNES_DIVERGED_FUNCTION_DOMAIN = -1, /* the new x location passed the function is not in the domain of F */
326184914b5SBarry Smith SNES_DIVERGED_FUNCTION_COUNT = -2,
32746a9e3ceSBarry Smith SNES_DIVERGED_LINEAR_SOLVE = -3, /* the linear solve failed */
32876c63389SBarry Smith SNES_DIVERGED_FUNCTION_NANORINF = -4,
3295664c22fSDavid Wells SNES_DIVERGED_FNORM_NAN_DEPRECATED = -4,
330184914b5SBarry Smith SNES_DIVERGED_MAX_IT = -5,
331647a2e1fSBarry Smith SNES_DIVERGED_LINE_SEARCH = -6, /* the line search failed */
3321e633543SBarry Smith SNES_DIVERGED_INNER = -7, /* inner solve failed */
33358c775ebSBarry Smith SNES_DIVERGED_LOCAL_MIN = -8, /* || J^T b || is small, implies converged to local minimum of F() */
334e37c518bSBarry Smith SNES_DIVERGED_DTOL = -9, /* || F || > divtol*||F_initial|| */
33507b62357SFande Kong SNES_DIVERGED_JACOBIAN_DOMAIN = -10, /* Jacobian calculation does not make sense */
3361c6b2ff8SBarry Smith SNES_DIVERGED_TR_DELTA = -11,
337c78072a7SPatrick Sanan SNES_CONVERGED_TR_DELTA_DEPRECATED = -11,
338379ef8d2SAlexander SNES_DIVERGED_USER = -12, /* The user has indicated divergence for an arbitrary reason */
33976c63389SBarry Smith SNES_DIVERGED_OBJECTIVE_DOMAIN = -13,
34076c63389SBarry Smith SNES_DIVERGED_OBJECTIVE_NANORINF = -14,
341e37c518bSBarry Smith
3429371c9d4SSatish Balay SNES_CONVERGED_ITERATING = 0
3439371c9d4SSatish Balay } SNESConvergedReason;
344014dd563SJed Brown PETSC_EXTERN const char *const *SNESConvergedReasons;
345184914b5SBarry Smith
346c838673bSBarry Smith /*MC
34776c63389SBarry Smith SNES_CONVERGED_FNORM_ABS - $||F|| \le abstol$
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
35576c63389SBarry Smith SNES_CONVERGED_FNORM_RELATIVE - $||F|| \le rtol*||F(x_0)||$ where $x_0$ is the initial guess
356c838673bSBarry Smith
357c838673bSBarry Smith Level: beginner
358c838673bSBarry Smith
3591cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
360c838673bSBarry Smith M*/
361c838673bSBarry Smith
362c838673bSBarry Smith /*MC
36376c63389SBarry Smith SNES_CONVERGED_SNORM_RELATIVE - The 2-norm of the last step $\le stol * ||x||$ where $x$ is the current
36476c63389SBarry Smith solution and $stol$ is the 4th argument to `SNESSetTolerances()`
365c838673bSBarry Smith
366af27ebaaSBarry Smith Options Database Key:
367ae9be289SBarry Smith -snes_stol <stol> - the step tolerance
368ae9be289SBarry Smith
369c838673bSBarry Smith Level: beginner
370c838673bSBarry Smith
3711cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
372c838673bSBarry Smith M*/
373c838673bSBarry Smith
374c838673bSBarry Smith /*MC
375c838673bSBarry Smith SNES_DIVERGED_FUNCTION_COUNT - The user provided function has been called more times then the final
37687497f52SBarry Smith argument to `SNESSetTolerances()`
377c838673bSBarry Smith
378c838673bSBarry Smith Level: beginner
379c838673bSBarry Smith
3801cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
381c838673bSBarry Smith M*/
382c838673bSBarry Smith
383c838673bSBarry Smith /*MC
38487497f52SBarry Smith SNES_DIVERGED_DTOL - The norm of the function has increased by a factor of divtol set with `SNESSetDivergenceTolerance()`
385e37c518bSBarry Smith
386e37c518bSBarry Smith Level: beginner
387e37c518bSBarry Smith
3881cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESSetDivergenceTolerance()`
389e37c518bSBarry Smith M*/
390e37c518bSBarry Smith
391e37c518bSBarry Smith /*MC
39276c63389SBarry Smith SNES_DIVERGED_FUNCTION_NANORINF - the 2-norm of the current function evaluation is not-a-number (NaN) or infinity, this
39376c63389SBarry Smith is usually caused by a division of 0 by 0, or infinity. See `SNESSetFunctionDomainError()`
39476c63389SBarry Smith
39576c63389SBarry Smith Level: beginner
39676c63389SBarry Smith
39776c63389SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
39876c63389SBarry Smith M*/
39976c63389SBarry Smith
40076c63389SBarry Smith /*MC
40176c63389SBarry Smith SNES_DIVERGED_FUNCTION_DOMAIN - the function provided with `SNESSetFunction()` called `SNESSetFunctionDomainError()` and
40276c63389SBarry Smith the solver could not recoverer.
40376c63389SBarry Smith
40476c63389SBarry Smith Level: beginner
40576c63389SBarry Smith
40676c63389SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
40776c63389SBarry Smith M*/
40876c63389SBarry Smith
40976c63389SBarry Smith /*MC
41076c63389SBarry Smith SNES_DIVERGED_OBJECTIVE_DOMAIN - the function provided with `SNESSetObjective()` called `SNESSetObjectiveDomainError()` and
41176c63389SBarry Smith the solver could not recoverer.
41276c63389SBarry Smith
41376c63389SBarry Smith Level: beginner
41476c63389SBarry Smith
41576c63389SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
41676c63389SBarry Smith M*/
41776c63389SBarry Smith
41876c63389SBarry Smith /*MC
41976c63389SBarry Smith SNES_DIVERGED_JACOBIAN_DOMAIN - the function provided with `SNESSetJacobian()` called `SNESSetJacobianDomainError()`
420c838673bSBarry Smith
421c838673bSBarry Smith Level: beginner
422c838673bSBarry Smith
4231cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
424c838673bSBarry Smith M*/
425c838673bSBarry Smith
426c838673bSBarry Smith /*MC
427c838673bSBarry Smith SNES_DIVERGED_MAX_IT - SNESSolve() has reached the maximum number of iterations requested
428c838673bSBarry Smith
429c838673bSBarry Smith Level: beginner
430c838673bSBarry Smith
4311cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
432c838673bSBarry Smith M*/
433c838673bSBarry Smith
434c838673bSBarry Smith /*MC
43587497f52SBarry Smith SNES_DIVERGED_LINE_SEARCH - The line search has failed. This only occurs for a `SNES` solvers that use a line search
436c838673bSBarry Smith
437c838673bSBarry Smith Level: beginner
438c838673bSBarry Smith
4391cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`, `SNESLineSearch`
440c838673bSBarry Smith M*/
441c838673bSBarry Smith
442c838673bSBarry Smith /*MC
44346a9e3ceSBarry Smith SNES_DIVERGED_LOCAL_MIN - the algorithm seems to have stagnated at a local minimum that is not zero.
44487497f52SBarry Smith See the manual page for `SNESConvergedReason` for more details
445c838673bSBarry Smith
446c838673bSBarry Smith Level: beginner
447c838673bSBarry Smith
4481cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
449c838673bSBarry Smith M*/
450c838673bSBarry Smith
451c838673bSBarry Smith /*MC
45287497f52SBarry Smith SNES_CONERGED_ITERATING - this only occurs if `SNESGetConvergedReason()` is called during the `SNESSolve()`
453c838673bSBarry Smith
454c838673bSBarry Smith Level: beginner
455c838673bSBarry Smith
4561cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `SNESConvergedReason`, `SNESSetTolerances()`
457c838673bSBarry Smith M*/
458c838673bSBarry Smith
45912651944SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetConvergenceTest(SNES, PetscErrorCode (*)(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *), void *, PetscCtxDestroyFn *);
4608d359177SBarry Smith PETSC_EXTERN PetscErrorCode SNESConvergedDefault(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
461e2a6519dSDmitry Karpeev PETSC_EXTERN PetscErrorCode SNESConvergedSkip(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
462649ef022SMatthew Knepley PETSC_EXTERN PetscErrorCode SNESConvergedCorrectPressure(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetConvergedReason(SNES, SNESConvergedReason *);
464c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode SNESGetConvergedReasonString(SNES, const char **);
46533866048SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetConvergedReason(SNES, SNESConvergedReason);
466ddbbdb52SLois Curfman McInnes
SNESSkipConverged(void)467edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "SNESConvergedSkip()", ) static inline void SNESSkipConverged(void)
468d71ae5a4SJacob Faibussowitsch { /* never called */
4699371c9d4SSatish Balay }
4708ea1b3e6SJed Brown #define SNESSkipConverged (SNESSkipConverged, SNESConvergedSkip)
47111f088b5SMatthew G Knepley
4729bcc50f1SBarry Smith /*S
4738434afd1SBarry Smith SNESInitialGuessFn - A prototype of a `SNES` compute initial guess function that would be passed to `SNESSetComputeInitialGuess()`
4749bcc50f1SBarry Smith
4759bcc50f1SBarry Smith Calling Sequence:
4769bcc50f1SBarry Smith + snes - `SNES` context
4779bcc50f1SBarry Smith . u - output vector to contain initial guess
4789bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4799bcc50f1SBarry Smith
4809bcc50f1SBarry Smith Level: beginner
4819bcc50f1SBarry Smith
4828434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetComputeInitialGuess()`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESFunctionFn`
4839bcc50f1SBarry Smith S*/
484*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESInitialGuessFn(SNES snes, Vec u, PetscCtx ctx);
4859bcc50f1SBarry Smith
4869bcc50f1SBarry Smith /*S
4878434afd1SBarry Smith SNESFunctionFn - A prototype of a `SNES` evaluation function that would be passed to `SNESSetFunction()`
4889bcc50f1SBarry Smith
4899bcc50f1SBarry Smith Calling Sequence:
4909bcc50f1SBarry Smith + snes - `SNES` context
4919bcc50f1SBarry Smith . u - input vector
4929bcc50f1SBarry Smith . F - function vector
4939bcc50f1SBarry Smith - ctx - [optional] user-defined function context
4949bcc50f1SBarry Smith
4959bcc50f1SBarry Smith Level: beginner
4969bcc50f1SBarry Smith
4978434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
4989bcc50f1SBarry Smith S*/
499*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESFunctionFn(SNES snes, Vec u, Vec F, PetscCtx ctx);
5009bcc50f1SBarry Smith
5019bcc50f1SBarry Smith /*S
5028434afd1SBarry Smith SNESObjectiveFn - A prototype of a `SNES` objective evaluation function that would be passed to `SNESSetObjective()`
5039bcc50f1SBarry Smith
5049bcc50f1SBarry Smith Calling Sequence:
5059bcc50f1SBarry Smith + snes - `SNES` context
5069bcc50f1SBarry Smith . u - input vector
5079bcc50f1SBarry Smith . o - output value
5089bcc50f1SBarry Smith - ctx - [optional] user-defined function context
5099bcc50f1SBarry Smith
5109bcc50f1SBarry Smith Level: beginner
5119bcc50f1SBarry Smith
5128434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`, `SNESNGSFn`
5139bcc50f1SBarry Smith S*/
514*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESObjectiveFn(SNES snes, Vec u, PetscReal *o, PetscCtx ctx);
5159bcc50f1SBarry Smith
5169bcc50f1SBarry Smith /*S
5178434afd1SBarry Smith SNESJacobianFn - A prototype of a `SNES` Jacobian evaluation function that would be passed to `SNESSetJacobian()`
5189bcc50f1SBarry Smith
5199bcc50f1SBarry Smith Calling Sequence:
5209bcc50f1SBarry Smith + snes - the `SNES` context obtained from `SNESCreate()`
5219bcc50f1SBarry Smith . u - input vector
5229bcc50f1SBarry Smith . Amat - (approximate) Jacobian matrix
5239bcc50f1SBarry Smith . Pmat - matrix used to construct the preconditioner, often the same as `Amat`
5249bcc50f1SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine
5259bcc50f1SBarry Smith
5269bcc50f1SBarry Smith Level: beginner
5279bcc50f1SBarry Smith
5288434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESNGSFn`
5299bcc50f1SBarry Smith S*/
530*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESJacobianFn(SNES snes, Vec u, Mat Amat, Mat Pmat, PetscCtx ctx);
5319bcc50f1SBarry Smith
5329bcc50f1SBarry Smith /*S
5338434afd1SBarry Smith SNESNGSFn - A prototype of a `SNES` nonlinear Gauss-Seidel function that would be passed to `SNESSetNGS()`
5349bcc50f1SBarry Smith
5359bcc50f1SBarry Smith Calling Sequence:
5369bcc50f1SBarry Smith + snes - the `SNES` context obtained from `SNESCreate()`
5379bcc50f1SBarry Smith . u - the current solution, updated in place
538dd8e379bSPierre Jolivet . b - the right-hand side vector (which may be `NULL`)
5399bcc50f1SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine
5409bcc50f1SBarry Smith
5419bcc50f1SBarry Smith Level: beginner
5429bcc50f1SBarry Smith
5438434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESGetJacobian()`, `SNESFunctionFn`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESJacobianFn`
5449bcc50f1SBarry Smith S*/
545*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESNGSFn(SNES snes, Vec u, Vec b, PetscCtx ctx);
5469bcc50f1SBarry Smith
54753e5d35bSStefano Zampini /*S
54853e5d35bSStefano Zampini SNESUpdateFn - A prototype of a `SNES` update function that would be passed to `SNESSetUpdate()`
54953e5d35bSStefano Zampini
55053e5d35bSStefano Zampini Calling Sequence:
55153e5d35bSStefano Zampini + snes - `SNES` context
55253e5d35bSStefano Zampini - step - the current iteration index
55353e5d35bSStefano Zampini
55453e5d35bSStefano Zampini Level: advanced
55553e5d35bSStefano Zampini
55653e5d35bSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESSetUpdate()`
55753e5d35bSStefano Zampini S*/
55821bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESUpdateFn(SNES snes, PetscInt step);
55953e5d35bSStefano Zampini
560b67197daSBarry Smith /* --------- Solving systems of nonlinear equations --------------- */
561*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetFunction(SNES, Vec, SNESFunctionFn *, PetscCtx);
562*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetFunction(SNES, Vec *, SNESFunctionFn **, PetscCtxRt);
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESComputeFunction(SNES, Vec, Vec);
564bbc1464cSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeMFFunction(SNES, Vec, Vec);
56525acbd8eSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESSetInitialFunction(SNES, Vec);
56632f3f7c2SPeter Brune
567*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetJacobian(SNES, Mat, Mat, SNESJacobianFn *, PetscCtx);
568*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetJacobian(SNES, Mat *, Mat *, SNESJacobianFn **, PetscCtxRt);
5698434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESObjectiveComputeFunctionDefaultFD;
5708434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefault;
5718434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESComputeJacobianDefaultColor;
57278e7fe0eSHong Zhang PETSC_EXTERN PetscErrorCode SNESPruneJacobianColor(SNES, Mat, Mat);
573*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetComputeInitialGuess(SNES, SNESInitialGuessFn *, PetscCtx);
574*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetPicard(SNES, Vec, SNESFunctionFn *, Mat, Mat, SNESJacobianFn *, PetscCtx);
575*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetPicard(SNES, Vec *, SNESFunctionFn **, Mat *, Mat *, SNESJacobianFn **, PetscCtxRt);
5768434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeFunction;
5778434afd1SBarry Smith PETSC_EXTERN SNESFunctionFn SNESPicardComputeMFFunction;
5788434afd1SBarry Smith PETSC_EXTERN SNESJacobianFn SNESPicardComputeJacobian;
57917bae607SBarry Smith
580*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetObjective(SNES, SNESObjectiveFn *, PetscCtx);
581*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetObjective(SNES, SNESObjectiveFn **, PetscCtxRt);
5822a4ee8f2SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeObjective(SNES, Vec, PetscReal *);
5832a4ee8f2SPeter Brune
58453e5d35bSStefano Zampini PETSC_EXTERN PetscErrorCode SNESSetUpdate(SNES, SNESUpdateFn *);
58553e5d35bSStefano Zampini
586534ebe21SPeter Brune /*E
58716a05f60SBarry Smith SNESNormSchedule - Frequency with which the norm is computed during a nonliner solve
588534ebe21SPeter Brune
589dc4c0fb0SBarry Smith Values:
590dc4c0fb0SBarry Smith + `SNES_NORM_DEFAULT` - use the default behavior for the current `SNESType`
591dc4c0fb0SBarry Smith . `SNES_NORM_NONE` - avoid all norm computations
592dc4c0fb0SBarry Smith . `SNES_NORM_ALWAYS` - compute the norms whenever possible
593dc4c0fb0SBarry Smith . `SNES_NORM_INITIAL_ONLY` - compute the norm only when the algorithm starts
594dc4c0fb0SBarry Smith . `SNES_NORM_FINAL_ONLY` - compute the norm only when the algorithm finishes
595dc4c0fb0SBarry Smith - `SNES_NORM_INITIAL_FINAL_ONLY` - compute the norm at the start and end of the algorithm
596534ebe21SPeter Brune
59716a05f60SBarry Smith Level: advanced
59816a05f60SBarry Smith
599534ebe21SPeter Brune Notes:
600dc4c0fb0SBarry Smith Support for these is highly dependent on the solver.
601dc4c0fb0SBarry Smith
602dc4c0fb0SBarry Smith Some options limit the convergence tests that can be used.
603dc4c0fb0SBarry Smith
604dc4c0fb0SBarry Smith The `SNES_NORM_NONE` option is most commonly used when the nonlinear solver is being used as a smoother, for example for `SNESFAS`
605dc4c0fb0SBarry Smith
606534ebe21SPeter Brune This is primarily used to turn off extra norm and function computation
607534ebe21SPeter Brune when the solvers are composed.
608534ebe21SPeter Brune
6091cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
610db781477SPatrick Sanan `KSPSetConvergenceTest()`, `KSPSetPCSide()`
611534ebe21SPeter Brune E*/
6129371c9d4SSatish Balay typedef enum {
6139371c9d4SSatish Balay SNES_NORM_DEFAULT = -1,
614fdacfa88SPeter Brune SNES_NORM_NONE = 0,
615365a6726SPeter Brune SNES_NORM_ALWAYS = 1,
616fdacfa88SPeter Brune SNES_NORM_INITIAL_ONLY = 2,
617fdacfa88SPeter Brune SNES_NORM_FINAL_ONLY = 3,
6189371c9d4SSatish Balay SNES_NORM_INITIAL_FINAL_ONLY = 4
6199371c9d4SSatish Balay } SNESNormSchedule;
620365a6726SPeter Brune PETSC_EXTERN const char *const *const SNESNormSchedules;
6211957e957SBarry Smith
622534ebe21SPeter Brune /*MC
62316a05f60SBarry Smith SNES_NORM_NONE - Don't compute function and its L2 norm when possible
624534ebe21SPeter Brune
625534ebe21SPeter Brune Level: advanced
626534ebe21SPeter Brune
627dc4c0fb0SBarry Smith Note:
628534ebe21SPeter Brune This is most useful for stationary solvers with a fixed number of iterations used as smoothers.
629534ebe21SPeter Brune
6301cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_DEFAULT`
631534ebe21SPeter Brune M*/
632534ebe21SPeter Brune
633fdacfa88SPeter Brune /*MC
634365a6726SPeter Brune SNES_NORM_ALWAYS - Compute the function and its L2 norm at each iteration.
635534ebe21SPeter Brune
636fdacfa88SPeter Brune Level: advanced
637fdacfa88SPeter Brune
638dc4c0fb0SBarry Smith Note:
639fdacfa88SPeter Brune Most solvers will use this no matter what norm type is passed to them.
640fdacfa88SPeter Brune
6411cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_NONE`
642fdacfa88SPeter Brune M*/
643534ebe21SPeter Brune
644534ebe21SPeter Brune /*MC
645534ebe21SPeter Brune SNES_NORM_INITIAL_ONLY - Compute the function and its L2 at iteration 0, but do not update it.
646534ebe21SPeter Brune
647534ebe21SPeter Brune Level: advanced
648534ebe21SPeter Brune
649534ebe21SPeter Brune Notes:
650dc4c0fb0SBarry Smith This method is useful in composed methods, when a true solution might actually be found before `SNESSolve()` is called.
651534ebe21SPeter Brune This option enables the solve to abort on the zeroth iteration if this is the case.
652534ebe21SPeter Brune
653534ebe21SPeter Brune For solvers that require the computation of the L2 norm of the function as part of the method, this merely cancels
654534ebe21SPeter Brune the norm computation at the last iteration (if possible).
655534ebe21SPeter Brune
6561cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_FINAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
657534ebe21SPeter Brune M*/
658534ebe21SPeter Brune
659534ebe21SPeter Brune /*MC
660534ebe21SPeter Brune SNES_NORM_FINAL_ONLY - Compute the function and its L2 norm on only the final iteration.
661534ebe21SPeter Brune
662534ebe21SPeter Brune Level: advanced
663534ebe21SPeter Brune
664dc4c0fb0SBarry Smith Note:
665534ebe21SPeter Brune For solvers that require the computation of the L2 norm of the function as part of the method, behaves
66687497f52SBarry Smith exactly as `SNES_NORM_DEFAULT`. This method is useful when the function is gotten after `SNESSolve()` and
667534ebe21SPeter Brune used in subsequent computation for methods that do not need the norm computed during the rest of the
668534ebe21SPeter Brune solution procedure.
669534ebe21SPeter Brune
6701cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_INITIAL_ONLY`, `SNES_NORM_INITIAL_FINAL_ONLY`
671534ebe21SPeter Brune M*/
672534ebe21SPeter Brune
673534ebe21SPeter Brune /*MC
674534ebe21SPeter Brune SNES_NORM_INITIAL_FINAL_ONLY - Compute the function and its L2 norm on only the initial and final iterations.
675534ebe21SPeter Brune
676534ebe21SPeter Brune Level: advanced
677534ebe21SPeter Brune
678dc4c0fb0SBarry Smith Note:
67987497f52SBarry Smith This method combines the benefits of `SNES_NORM_INITIAL_ONLY` and `SNES_NORM_FINAL_ONLY`.
680534ebe21SPeter Brune
6811cc06b55SBarry Smith .seealso: [](ch_snes), `SNESNormSchedule`, `SNES`, `SNESSetNormSchedule()`, `SNES_NORM_SNES_NORM_INITIAL_ONLY`, `SNES_NORM_FINAL_ONLY`
682534ebe21SPeter Brune M*/
683534ebe21SPeter Brune
684365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetNormSchedule(SNES, SNESNormSchedule);
685365a6726SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetNormSchedule(SNES, SNESNormSchedule *);
686c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESSetFunctionNorm(SNES, PetscReal);
687c5ce4427SMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESGetFunctionNorm(SNES, PetscReal *);
688c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetUpdateNorm(SNES, PetscReal *);
689c1e67a49SFande Kong PETSC_EXTERN PetscErrorCode SNESGetSolutionNorm(SNES, PetscReal *);
690534ebe21SPeter Brune
69147073ea2SPeter Brune /*E
69247073ea2SPeter Brune SNESFunctionType - Type of function computed
69347073ea2SPeter Brune
694dc4c0fb0SBarry Smith Values:
695dc4c0fb0SBarry Smith + `SNES_FUNCTION_DEFAULT` - the default behavior for the current `SNESType`
696dc4c0fb0SBarry Smith . `SNES_FUNCTION_UNPRECONDITIONED` - the original function provided
697dc4c0fb0SBarry Smith - `SNES_FUNCTION_PRECONDITIONED` - the modification of the function by the preconditioner
69847073ea2SPeter Brune
69916a05f60SBarry Smith Level: advanced
70016a05f60SBarry Smith
701dc4c0fb0SBarry Smith Note:
702dc4c0fb0SBarry Smith Support for these is dependent on the solver.
703dc4c0fb0SBarry Smith
7041cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPSetNormType()`,
705db781477SPatrick Sanan `KSPSetConvergenceTest()`, `KSPSetPCSide()`
70647073ea2SPeter Brune E*/
7079371c9d4SSatish Balay typedef enum {
7089371c9d4SSatish Balay SNES_FUNCTION_DEFAULT = -1,
70947073ea2SPeter Brune SNES_FUNCTION_UNPRECONDITIONED = 0,
7109371c9d4SSatish Balay SNES_FUNCTION_PRECONDITIONED = 1
7119371c9d4SSatish Balay } SNESFunctionType;
71247073ea2SPeter Brune PETSC_EXTERN const char *const *const SNESFunctionTypes;
71347073ea2SPeter Brune
71447073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESSetFunctionType(SNES, SNESFunctionType);
71547073ea2SPeter Brune PETSC_EXTERN PetscErrorCode SNESGetFunctionType(SNES, SNESFunctionType *);
716f1c6b773SPeter Brune
717*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNGS(SNES, SNESNGSFn *, PetscCtx);
718*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNGS(SNES, SNESNGSFn **, PetscCtxRt);
719be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeNGS(SNES, Vec, Vec);
720b6266c6eSPeter Brune
721be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetSweeps(SNES, PetscInt);
722be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetSweeps(SNES, PetscInt *);
723be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSSetTolerances(SNES, PetscReal, PetscReal, PetscReal, PetscInt);
724be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESNGSGetTolerances(SNES, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
725badc63e7SPeter Brune
7264fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESSetAlwaysComputesFinalResidual(SNES, PetscBool);
7274fc747eaSLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESGetAlwaysComputesFinalResidual(SNES, PetscBool *);
7284fc747eaSLawrence Mitchell
729*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESShellGetContext(SNES, PetscCtxRt);
730*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESShellSetContext(SNES, PetscCtx);
731014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESShellSetSolve(SNES, PetscErrorCode (*)(SNES, Vec));
732646217ecSPeter Brune
733c5ae4b9aSBarry Smith /* --------- Routines specifically for line search methods --------------- */
734c5ae4b9aSBarry Smith
735872b6db9SPeter Brune /*S
73687497f52SBarry Smith SNESLineSearch - Abstract PETSc object that manages line-search operations for nonlinear solvers
7379e764e56SPeter Brune
7389e764e56SPeter Brune Level: beginner
7399e764e56SPeter Brune
7401cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNES`
7419e764e56SPeter Brune S*/
742907376e6SBarry Smith typedef struct _p_LineSearch *SNESLineSearch;
7439e764e56SPeter Brune
7449e764e56SPeter Brune /*J
7450b4b7b1cSBarry Smith SNESLineSearchType - String with the name of a PETSc line search method `SNESLineSearch`. Provides all the linesearches for the nonlinear solvers, `SNES`,
7460b4b7b1cSBarry Smith in PETSc.
7479e764e56SPeter Brune
748ceaaa498SBarry Smith Values:
749ceaaa498SBarry Smith + `SNESLINESEARCHBASIC` - (or equivalently `SNESLINESEARCHNONE`) Simple damping line search, defaults to using the full Newton step
750a99ef635SJonas Heinzmann . `SNESLINESEARCHBT` - Backtracking line search over the L2 norm of the function or an objective function
751a99ef635SJonas Heinzmann . `SNESLINESEARCHSECANT` - Secant line search over the L2 norm of the function or an objective function
7520b4b7b1cSBarry Smith . `SNESLINESEARCHCP` - Critical point secant line search assuming $F(x) = \nabla G(x)$ for some unknown $G(x)$
753ceaaa498SBarry Smith . `SNESLINESEARCHNLEQERR` - Affine-covariant error-oriented linesearch
754a99ef635SJonas Heinzmann - `SNESLINESEARCHBISECTION` - bisection line search for a root in the directional derivative
755ceaaa498SBarry Smith - `SNESLINESEARCHSHELL` - User provided `SNESLineSearch` implementation
756ceaaa498SBarry Smith
75795bd0b28SBarry Smith Level: beginner
75895bd0b28SBarry Smith
7590b4b7b1cSBarry Smith Note:
7600b4b7b1cSBarry Smith Use `SNESLineSearchSetType()` or the options database key `-snes_linesearch_type` to set
7610b4b7b1cSBarry Smith the specific line search algorithm to use with a given `SNES` object. Not all `SNESType` can utilize a line search.
7620b4b7b1cSBarry Smith
7631cc06b55SBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchSetType()`, `SNES`
7649e764e56SPeter Brune J*/
76519fd82e9SBarry Smith typedef const char *SNESLineSearchType;
7669e764e56SPeter Brune #define SNESLINESEARCHBT "bt"
767d4c6564cSPatrick Farrell #define SNESLINESEARCHNLEQERR "nleqerr"
768c87759e9SPeter Brune #define SNESLINESEARCHBASIC "basic"
7690b00b554SBarry Smith #define SNESLINESEARCHNONE "none"
770a99ef635SJonas Heinzmann #define SNESLINESEARCHSECANT "secant"
771392273beSMatthew G. Knepley #define SNESLINESEARCHL2 PETSC_DEPRECATED_MACRO(3, 24, 0, "SNESLINESEARCHSECANT", ) "secant"
7729e764e56SPeter Brune #define SNESLINESEARCHCP "cp"
7739e764e56SPeter Brune #define SNESLINESEARCHSHELL "shell"
774b5badacbSBarry Smith #define SNESLINESEARCHNCGLINEAR "ncglinear"
775b9d635d7SJonas Heinzmann #define SNESLINESEARCHBISECTION "bisection"
7769e764e56SPeter Brune
777140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESList;
778014dd563SJed Brown PETSC_EXTERN PetscClassId SNESLINESEARCH_CLASSID;
779140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList SNESLineSearchList;
7809e764e56SPeter Brune
781b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_LINEAR 1
782b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_QUADRATIC 2
783b000cd8dSPeter Brune #define SNES_LINESEARCH_ORDER_CUBIC 3
7849e764e56SPeter Brune
7859bcc50f1SBarry Smith /*S
7868434afd1SBarry Smith SNESLineSearchVIProjectFn - A prototype of a `SNES` function that projects a vector onto the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
7879bcc50f1SBarry Smith
7889bcc50f1SBarry Smith Calling Sequence:
7899bcc50f1SBarry Smith + snes - `SNES` context
7909bcc50f1SBarry Smith - u - the vector to project to the bounds
7919bcc50f1SBarry Smith
7929bcc50f1SBarry Smith Level: advanced
7939bcc50f1SBarry Smith
7949bcc50f1SBarry Smith Note:
7958434afd1SBarry Smith The deprecated `SNESLineSearchVIProjectFunc` still works as a replacement for `SNESLineSearchVIProjectFn` *.
7969bcc50f1SBarry Smith
7979bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
7989bcc50f1SBarry Smith S*/
79921bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchVIProjectFn(SNES snes, Vec u);
80021bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVIProjectFn *SNESLineSearchVIProjectFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVIProjectFn*", );
8019bcc50f1SBarry Smith
8029bcc50f1SBarry Smith /*S
8038434afd1SBarry Smith SNESLineSearchVIProjectFn - A prototype of a `SNES` function that computes the norm of the active set variables in a vector in a VI solve,
8049bcc50f1SBarry Smith passed to `SNESLineSearchSetVIFunctions()`
8059bcc50f1SBarry Smith
8069bcc50f1SBarry Smith Calling Sequence:
8079bcc50f1SBarry Smith + snes - `SNES` context
8089bcc50f1SBarry Smith . f - the vector to compute the norm of
8099bcc50f1SBarry Smith . u - the current solution, entries that are on the VI bounds are ignored
8109bcc50f1SBarry Smith - fnorm - the resulting norm
8119bcc50f1SBarry Smith
8129bcc50f1SBarry Smith Level: advanced
8139bcc50f1SBarry Smith
8149bcc50f1SBarry Smith Note:
8158434afd1SBarry Smith The deprecated `SNESLineSearchVINormFunc` still works as a replacement for `SNESLineSearchVINormFn` *.
8169bcc50f1SBarry Smith
8179bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch`
8189bcc50f1SBarry Smith S*/
81921bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchVINormFn(SNES snes, Vec f, Vec u, PetscReal *fnorm);
82021bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchVINormFn *SNESLineSearchVINormFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchVINormFnn*", );
8219bcc50f1SBarry Smith
822d5def619SJonas Heinzmann /*S
823d5def619SJonas Heinzmann SNESLineSearchVIDirDerivFn - A prototype of a `SNES` function that computes the directional derivative considering the VI bounds, passed to `SNESLineSearchSetVIFunctions()`
824d5def619SJonas Heinzmann
825d5def619SJonas Heinzmann Calling Sequence:
826d5def619SJonas Heinzmann + snes - `SNES` context
827d5def619SJonas Heinzmann . f - the function vector to compute the directional derivative with
828d5def619SJonas Heinzmann . u - the current solution, entries that are on the VI bounds are ignored
829d5def619SJonas Heinzmann . y - the direction to compute the directional derivative
830d5def619SJonas Heinzmann - fty - the resulting directional derivative
831d5def619SJonas Heinzmann
832d5def619SJonas Heinzmann Level: advanced
833d5def619SJonas Heinzmann
834d5def619SJonas Heinzmann .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchVIProjectFn`, `SNESLineSearchVIProjectFn`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetVIFunctions()`
835d5def619SJonas Heinzmann S*/
83621bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchVIDirDerivFn(SNES snes, Vec f, Vec u, Vec y, PetscScalar *fty);
837d5def619SJonas Heinzmann
83821bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchApplyFn(SNESLineSearch);
83921bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchApplyFn *SNESLineSearchApplyFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
840*2a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode SNESLineSearchShellApplyFn(SNESLineSearch, PetscCtx);
84121bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef SNESLineSearchShellApplyFn *SNESLineSearchUserFunc PETSC_DEPRECATED_TYPEDEF(3, 21, 0, "SNESLineSearchApplyFn*", );
8429e764e56SPeter Brune
843014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate(MPI_Comm, SNESLineSearch *);
844014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchReset(SNESLineSearch);
845014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchView(SNESLineSearch, PetscViewer);
846014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *);
847a80ff896SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetType(SNESLineSearch, SNESLineSearchType *);
84819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetType(SNESLineSearch, SNESLineSearchType);
849014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch);
850ed07d7d7SPeter Brune PETSC_EXTERN PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch, PetscErrorCode (*)(SNES, Vec, Vec));
851014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetUp(SNESLineSearch);
852014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchApply(SNESLineSearch, Vec, Vec, PetscReal *, Vec);
853014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch, Vec, Vec, PetscBool *);
854014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *);
855fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch, PetscInt);
8569e764e56SPeter Brune
8579bd66eb0SPeter Brune /* set the functions for precheck and postcheck */
85886d74e61SPeter Brune
859*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, PetscBool *, PetscCtx), PetscCtx);
860*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, PetscCtx), PetscCtx);
8619e764e56SPeter Brune
862*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, PetscBool *, PetscCtx), PetscCtxRt);
863*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch, PetscErrorCode (**)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, PetscCtx), PetscCtxRt);
8649e764e56SPeter Brune
8659bd66eb0SPeter Brune /* set the functions for VI-specific line search operations */
8669bd66eb0SPeter Brune
867d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn *, SNESLineSearchVINormFn *, SNESLineSearchVIDirDerivFn *);
868d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch, SNESLineSearchVIProjectFn **, SNESLineSearchVINormFn **, SNESLineSearchVIDirDerivFn **);
8699bd66eb0SPeter Brune
8709e764e56SPeter Brune /* pointers to the associated SNES in order to be able to get the function evaluation out */
871014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch, SNES);
872014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch, SNES *);
8739e764e56SPeter Brune
8749e764e56SPeter Brune /* set and get the parameters and vectors */
875014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscInt *);
876014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscInt);
8779e764e56SPeter Brune
878*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch, Vec, Vec, PetscBool *, PetscCtx);
87986d74e61SPeter Brune
880014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch, PetscReal *);
881014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch, PetscReal);
8829e764e56SPeter Brune
883014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch, PetscReal *);
884014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch, PetscReal);
8859e764e56SPeter Brune
886ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch, PetscInt *);
887ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch, PetscInt);
88859405d5eSPeter Brune
889422a814eSBarry Smith /*E
890dc4c0fb0SBarry Smith SNESLineSearchReason - indication if the line search has succeeded or failed and why
891422a814eSBarry Smith
892dc4c0fb0SBarry Smith Values:
893dc4c0fb0SBarry Smith + `SNES_LINESEARCH_SUCCEEDED` - the line search succeeded
894dc4c0fb0SBarry Smith . `SNES_LINESEARCH_FAILED_NANORINF` - a not a number of infinity appeared in the computions
89576c63389SBarry Smith . `SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN` - the function was evaluated outside of its domain, see `SNESSetFunctionDomainError()`
89676c63389SBarry Smith . `SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN`- the objective function was evaluated outside of its domain, see `SNESSetObjectiveDomainError()`
89776c63389SBarry Smith . `SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN` - the Jacobian was evaluated outside of its domain, see `SNESSetJacobianDomainError()`
898dc4c0fb0SBarry Smith . `SNES_LINESEARCH_FAILED_REDUCT` - the linear search failed to get the requested decrease in its norm or objective
899dc4c0fb0SBarry Smith . `SNES_LINESEARCH_FAILED_USER` - used by `SNESLINESEARCHNLEQERR` to indicate the user changed the search direction inappropriately
900dc4c0fb0SBarry Smith - `SNES_LINESEARCH_FAILED_FUNCTION` - indicates the maximum number of function evaluations allowed has been surpassed, `SNESConvergedReason` is also
901dc4c0fb0SBarry Smith set to `SNES_DIVERGED_FUNCTION_COUNT`
902422a814eSBarry Smith
90316a05f60SBarry Smith Level: intermediate
90416a05f60SBarry Smith
905dc4c0fb0SBarry Smith Developer Note:
90676c63389SBarry Smith Some of these reasons overlap with values of `SNESConvergedReason`. It is possibly a better design to have `SNESConvergedReaon` alone used also for indicating line
90776c63389SBarry Smith search failures.
908422a814eSBarry Smith
9091cc06b55SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSolve()`, `SNESGetConvergedReason()`, `KSPConvergedReason`, `SNESSetConvergenceTest()`,
910dc4c0fb0SBarry Smith `SNESSetFunctionDomainError()` and `SNESSetJacobianDomainError()`
911422a814eSBarry Smith E*/
9129371c9d4SSatish Balay typedef enum {
9139371c9d4SSatish Balay SNES_LINESEARCH_SUCCEEDED,
914422a814eSBarry Smith SNES_LINESEARCH_FAILED_NANORINF,
91576c63389SBarry Smith SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN,
91676c63389SBarry Smith SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN,
91776c63389SBarry Smith SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN,
918d5b43468SJose E. Roman SNES_LINESEARCH_FAILED_REDUCT, /* INSUFFICIENT REDUCTION */
919e9b602ebSSatish Balay SNES_LINESEARCH_FAILED_USER,
9209371c9d4SSatish Balay SNES_LINESEARCH_FAILED_FUNCTION
9219371c9d4SSatish Balay } SNESLineSearchReason;
922422a814eSBarry Smith
923422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetReason(SNESLineSearch, SNESLineSearchReason *);
924422a814eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetReason(SNESLineSearch, SNESLineSearchReason);
9259e764e56SPeter Brune
926014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch, Vec *, Vec *, Vec *, Vec *, Vec *);
927014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch, Vec, Vec, Vec, Vec, Vec);
9289e764e56SPeter Brune
929014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch, PetscReal *, PetscReal *, PetscReal *);
930014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch, PetscReal, PetscReal, PetscReal);
931014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch);
932014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch, PetscBool);
9339e764e56SPeter Brune
934dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitor(SNESLineSearch);
935*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch, PetscErrorCode (*)(SNESLineSearch, PetscCtx), PetscCtx, PetscCtxDestroyFn *);
936d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch, const char[], const char[], const char[], PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*)(SNESLineSearch, PetscViewerAndFormat *));
937dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch);
938dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch, PetscViewer);
939dcf2fd19SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch, PetscViewer *);
940d12e167eSBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch, PetscViewerAndFormat *);
9419e764e56SPeter Brune
942ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch, const char[]);
943ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch, const char *[]);
9449e764e56SPeter Brune
9459e764e56SPeter Brune /* Shell interface functions */
946*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellSetApply(SNESLineSearch, SNESLineSearchShellApplyFn *, PetscCtx);
947*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchShellGetApply(SNESLineSearch, SNESLineSearchShellApplyFn **, PetscCtxRt);
9489bcc50f1SBarry Smith
SNESLineSearchShellSetUserFunc(SNESLineSearch ls,SNESLineSearchShellApplyFn * f,PetscCtx ctx)949*2a8381b2SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellSetApply()", ) static inline PetscErrorCode SNESLineSearchShellSetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn *f, PetscCtx ctx)
9509bcc50f1SBarry Smith {
9519bcc50f1SBarry Smith return SNESLineSearchShellSetApply(ls, f, ctx);
9529bcc50f1SBarry Smith }
9539bcc50f1SBarry Smith
SNESLineSearchShellGetUserFunc(SNESLineSearch ls,SNESLineSearchShellApplyFn ** f,PetscCtxRt ctx)954*2a8381b2SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 21, 0, "SNESLinesearchShellGetApply()", ) static inline PetscErrorCode SNESLineSearchShellGetUserFunc(SNESLineSearch ls, SNESLineSearchShellApplyFn **f, PetscCtxRt ctx)
9559bcc50f1SBarry Smith {
9569bcc50f1SBarry Smith return SNESLineSearchShellGetApply(ls, f, ctx);
9579bcc50f1SBarry Smith }
9589e764e56SPeter Brune
9592f4102e2SPeter Brune /* BT interface functions */
960014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch, PetscReal);
961014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch, PetscReal *);
9622f4102e2SPeter Brune
9639e764e56SPeter Brune /*register line search types */
964bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode SNESLineSearchRegister(const char[], PetscErrorCode (*)(SNESLineSearch));
9659e764e56SPeter Brune
966720c9a41SShri Abhyankar /* Routines for VI solver */
967014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetVariableBounds(SNES, Vec, Vec);
968cf836535SBlaise Bourdin PETSC_EXTERN PetscErrorCode SNESVIGetVariableBounds(SNES, Vec *, Vec *);
969014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVISetComputeVariableBounds(SNES, PetscErrorCode (*)(SNES, Vec, Vec));
970014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetInactiveSet(SNES, IS *);
971014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIGetActiveSetIS(SNES, Vec, Vec, IS *);
972014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFnorm(SNES, Vec, Vec, PetscReal *);
973d5def619SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESVIComputeInactiveSetFtY(SNES, Vec, Vec, Vec, PetscScalar *);
974*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESVISetRedundancyCheck(SNES, PetscErrorCode (*)(SNES, IS, IS *, PetscCtx), PetscCtx);
975f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeMeritFunction(Vec, PetscReal *, PetscReal *);
976f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode SNESVIComputeFunction(SNES, Vec, Vec, void *);
977f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMSetVI(DM, IS);
978f6dfbefdSBarry Smith PETSC_EXTERN PetscErrorCode DMDestroyVI(DM);
979f5ea5bd2SShri Abhyankar
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTestLocalMin(SNES);
981da9b6338SBarry Smith
982eef9c623SLois Curfman McInnes /* Should this routine be private? */
983d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeJacobian(SNES, Vec, Mat, Mat);
984cbf8f02cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode SNESTestJacobian(SNES, PetscReal *, PetscReal *);
985494a190aSStefano Zampini PETSC_EXTERN PetscErrorCode SNESTestFunction(SNES);
986ddbbdb52SLois Curfman McInnes
987014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetDM(SNES, DM);
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESGetDM(SNES, DM *);
989be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPC(SNES, SNES);
990be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPC(SNES, SNES *);
9913ad1a0b9SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESHasNPC(SNES, PetscBool *);
992be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESApplyNPC(SNES, Vec, Vec, Vec);
993be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCFunction(SNES, Vec, PetscReal *);
994be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESComputeFunctionDefaultNPC(SNES, Vec, Vec);
995be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESSetNPCSide(SNES, PCSide);
996be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESGetNPCSide(SNES, PCSide *);
9977601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESSetLineSearch(SNES, SNESLineSearch);
9987601faf0SJed Brown PETSC_EXTERN PetscErrorCode SNESGetLineSearch(SNES, SNESLineSearch *);
9996c699258SBarry Smith
SNESGetSNESLineSearch(SNES snes,SNESLineSearch * ls)1000edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESGetLineSearch()", ) static inline PetscErrorCode SNESGetSNESLineSearch(SNES snes, SNESLineSearch *ls)
1001d71ae5a4SJacob Faibussowitsch {
10029371c9d4SSatish Balay return SNESGetLineSearch(snes, ls);
10039371c9d4SSatish Balay }
SNESSetSNESLineSearch(SNES snes,SNESLineSearch ls)1004edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 4, 0, "SNESSetLineSearch()", ) static inline PetscErrorCode SNESSetSNESLineSearch(SNES snes, SNESLineSearch ls)
1005d71ae5a4SJacob Faibussowitsch {
10069371c9d4SSatish Balay return SNESSetLineSearch(snes, ls);
10079371c9d4SSatish Balay }
10088b7b3213SJed Brown
1009014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESSetUpMatrices(SNES);
1010*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunction(DM, SNESFunctionFn *, PetscCtx);
1011*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetFunction(DM, SNESFunctionFn **, PetscCtxRt);
101249abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionContextDestroy(DM, PetscCtxDestroyFn *);
1013*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetMFFunction(DM, SNESFunctionFn *, PetscCtx);
1014*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetNGS(DM, SNESNGSFn *, PetscCtx);
1015*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetNGS(DM, SNESNGSFn **, PetscCtxRt);
1016*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobian(DM, SNESJacobianFn *, PetscCtx);
1017*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetJacobian(DM, SNESJacobianFn **, PetscCtxRt);
101849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianContextDestroy(DM, PetscCtxDestroyFn *);
1019*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetPicard(DM, SNESFunctionFn *, SNESJacobianFn *, PetscCtx);
1020*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetPicard(DM, SNESFunctionFn **, SNESJacobianFn **, PetscCtxRt);
1021*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetObjective(DM, SNESObjectiveFn *, PetscCtx);
1022*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetObjective(DM, SNESObjectiveFn **, PetscCtxRt);
10234dd50a75SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMSNES(DM, DM);
10246cab3a1bSJed Brown
102521bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionFn(DMDALocalInfo *, void *, void *, void *);
102621bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianFn(DMDALocalInfo *, void *, Mat, Mat, void *);
102721bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveFn(DMDALocalInfo *, void *, PetscReal *, void *);
1028c9d099b5SPeter Brune
102921bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESFunctionVecFn(DMDALocalInfo *, Vec, Vec, void *);
103021bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESJacobianVecFn(DMDALocalInfo *, Vec, Mat, Mat, void *);
103121bb954dSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode DMDASNESObjectiveVecFn(DMDALocalInfo *, Vec, PetscReal *, void *);
10325505f7afSJunchao Zhang
10338434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocal(DM, InsertMode, DMDASNESFunctionFn *, void *);
10348434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocal(DM, DMDASNESJacobianFn *, void *);
10358434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocal(DM, DMDASNESObjectiveFn *, void *);
10368434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetPicardLocal(DM, InsertMode, DMDASNESFunctionFn *, DMDASNESJacobianFn, void *);
10376cab3a1bSJed Brown
10388434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetFunctionLocalVec(DM, InsertMode, DMDASNESFunctionVecFn *, void *);
10398434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetJacobianLocalVec(DM, DMDASNESJacobianVecFn *, void *);
10408434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDASNESSetObjectiveLocalVec(DM, DMDASNESObjectiveVecFn *, void *);
10415505f7afSJunchao Zhang
1042*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetBoundaryLocal(DM, PetscErrorCode (*)(DM, Vec, PetscCtx), PetscCtx);
1043*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetObjectiveLocal(DM, PetscErrorCode (*)(DM, Vec, PetscReal *, PetscCtx), PetscCtx);
1044*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetFunctionLocal(DM, PetscErrorCode (*)(DM, Vec, Vec, PetscCtx), PetscCtx);
1045*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESSetJacobianLocal(DM, PetscErrorCode (*)(DM, Vec, Mat, Mat, PetscCtx), PetscCtx);
1046*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetBoundaryLocal(DM, PetscErrorCode (**)(DM, Vec, PetscCtx), PetscCtxRt);
1047*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetObjectiveLocal(DM, PetscErrorCode (**)(DM, Vec, PetscReal *, PetscCtx), PetscCtxRt);
1048*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetFunctionLocal(DM, PetscErrorCode (**)(DM, Vec, Vec, PetscCtx), PetscCtxRt);
1049*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESGetJacobianLocal(DM, PetscErrorCode (**)(DM, Vec, Mat, Mat, PetscCtx), PetscCtxRt);
1050ff35dfedSBarry Smith
1051b665c654SMatthew G Knepley /* Routines for Multiblock solver */
1052014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetFields(SNES, const char[], PetscInt, const PetscInt *);
1053014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetIS(SNES, const char[], IS);
1054014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetBlockSize(SNES, PetscInt);
1055014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMultiblockSetType(SNES, PCCompositeType);
10564bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMultiblockGetSubSNES(SNES, PetscInt *, SNES *[]);
1057aeed3662SMatthew G Knepley
105837e1895aSJed Brown /*J
105987497f52SBarry Smith SNESMSType - String with the name of a PETSc `SNESMS` method.
106037e1895aSJed Brown
106137e1895aSJed Brown Level: intermediate
106237e1895aSJed Brown
10631cc06b55SBarry Smith .seealso: [](ch_snes), `SNESMS`, `SNESMSGetType()`, `SNESMSSetType()`, `SNES`
106437e1895aSJed Brown J*/
106519fd82e9SBarry Smith typedef const char *SNESMSType;
106637e1895aSJed Brown #define SNESMSM62 "m62"
106737e1895aSJed Brown #define SNESMSEULER "euler"
1068a97cb6bcSJed Brown #define SNESMSJAMESON83 "jameson83"
10693847c725SLisandro Dalcin #define SNESMSVLTP11 "vltp11"
1070a97cb6bcSJed Brown #define SNESMSVLTP21 "vltp21"
1071a97cb6bcSJed Brown #define SNESMSVLTP31 "vltp31"
1072a97cb6bcSJed Brown #define SNESMSVLTP41 "vltp41"
1073a97cb6bcSJed Brown #define SNESMSVLTP51 "vltp51"
1074a97cb6bcSJed Brown #define SNESMSVLTP61 "vltp61"
107537e1895aSJed Brown
107619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSRegister(SNESMSType, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[]);
10774bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESMSRegisterAll(void);
107857715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetType(SNES, SNESMSType *);
107919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSSetType(SNES, SNESMSType);
108057715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSGetDamping(SNES, PetscReal *);
108157715debSLisandro Dalcin PETSC_EXTERN PetscErrorCode SNESMSSetDamping(SNES, PetscReal);
1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSFinalizePackage(void);
1083607a6623SBarry Smith PETSC_EXTERN PetscErrorCode SNESMSInitializePackage(void);
1084014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESMSRegisterDestroy(void);
108537e1895aSJed Brown
1086ceaaa498SBarry Smith /*MC
1087ceaaa498SBarry Smith SNESNGMRESRestartType - the restart approach used by `SNESNGMRES`
108813a62661SPeter Brune
1089ceaaa498SBarry Smith Values:
1090ceaaa498SBarry Smith + `SNES_NGMRES_RESTART_NONE` - never restart
1091ceaaa498SBarry Smith . `SNES_NGMRES_RESTART_DIFFERENCE` - restart based upon difference criteria
1092ceaaa498SBarry Smith - `SNES_NGMRES_RESTART_PERIODIC` - restart after a fixed number of iterations
1093ceaaa498SBarry Smith
1094ceaaa498SBarry Smith Options Database Keys:
1095ceaaa498SBarry Smith + -snes_ngmres_restart_type <difference,periodic,none> - set the restart type
1096af27ebaaSBarry Smith - -snes_ngmres_restart <30> - sets the number of iterations before restart for periodic
1097ceaaa498SBarry Smith
109895bd0b28SBarry Smith Level: intermediate
109995bd0b28SBarry Smith
1100ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1101ceaaa498SBarry Smith `SNESNGMRESGetRestartType()`, `SNESNGMRESSelectType`
1102ceaaa498SBarry Smith M*/
110313a62661SPeter Brune typedef enum {
110413a62661SPeter Brune SNES_NGMRES_RESTART_NONE = 0,
110513a62661SPeter Brune SNES_NGMRES_RESTART_PERIODIC = 1,
11069371c9d4SSatish Balay SNES_NGMRES_RESTART_DIFFERENCE = 2
11079371c9d4SSatish Balay } SNESNGMRESRestartType;
11086a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESRestartTypes[];
110913a62661SPeter Brune
1110ceaaa498SBarry Smith /*MC
1111ceaaa498SBarry Smith SNESNGMRESSelectType - the approach used by `SNESNGMRES` to determine how the candidate solution and
1112ceaaa498SBarry Smith combined solution are used to create the next iterate.
1113ceaaa498SBarry Smith
1114ceaaa498SBarry Smith Values:
1115ceaaa498SBarry Smith + `SNES_NGMRES_SELECT_NONE` - choose the combined solution all the time
1116ceaaa498SBarry Smith . `SNES_NGMRES_SELECT_DIFFERENCE` - choose based upon the selection criteria
1117ceaaa498SBarry Smith - `SNES_NGMRES_SELECT_LINESEARCH` - choose based upon line search combination
1118ceaaa498SBarry Smith
1119ceaaa498SBarry Smith Options Database Key:
1120ceaaa498SBarry Smith . -snes_ngmres_select_type<difference,none,linesearch> - select type
1121ceaaa498SBarry Smith
112295bd0b28SBarry Smith Level: intermediate
112395bd0b28SBarry Smith
1124ceaaa498SBarry Smith .seealso: `SNES, `SNESNGMRES`, `SNESNGMRESSetSelectType()`, `SNESNGMRESGetSelectType()`, `SNESNGMRESSetRestartType()`,
1125ceaaa498SBarry Smith `SNESNGMRESGetRestartType()`, `SNESNGMRESRestartType`
1126ceaaa498SBarry Smith M*/
112713a62661SPeter Brune typedef enum {
112813a62661SPeter Brune SNES_NGMRES_SELECT_NONE = 0,
112913a62661SPeter Brune SNES_NGMRES_SELECT_DIFFERENCE = 1,
11309371c9d4SSatish Balay SNES_NGMRES_SELECT_LINESEARCH = 2
11319371c9d4SSatish Balay } SNESNGMRESSelectType;
11326a6fc655SJed Brown PETSC_EXTERN const char *const SNESNGMRESSelectTypes[];
113313a62661SPeter Brune
1134014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartType(SNES, SNESNGMRESRestartType);
1135014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNGMRESSetSelectType(SNES, SNESNGMRESSelectType);
113623b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESSetRestartFmRise(SNES, PetscBool);
113723b3e82cSAsbjørn Nilsen Riseth PETSC_EXTERN PetscErrorCode SNESNGMRESGetRestartFmRise(SNES, PetscBool *);
113813a62661SPeter Brune
1139ceaaa498SBarry Smith /*MC
1140ceaaa498SBarry Smith SNESNCGType - the conjugate update approach for `SNESNCG`
11410a844d1aSPeter Brune
1142ceaaa498SBarry Smith Values:
1143ceaaa498SBarry Smith + `SNES_NCG_FR` - Fletcher-Reeves update
1144ceaaa498SBarry Smith . `SNES_NCG_PRP` - Polak-Ribiere-Polyak update, the default and the only one that tolerates generalized search directions
1145ceaaa498SBarry Smith . `SNES_NCG_HS` - Hestenes-Steifel update
1146ceaaa498SBarry Smith . `SNES_NCG_DY` - Dai-Yuan update
1147ceaaa498SBarry Smith - `SNES_NCG_CD` - Conjugate Descent update
1148ceaaa498SBarry Smith
1149ceaaa498SBarry Smith Options Database Key:
1150ceaaa498SBarry Smith . -snes_ncg_type<fr,prp,hs,dy,cd> - select type
1151ceaaa498SBarry Smith
115295bd0b28SBarry Smith Level: intermediate
115395bd0b28SBarry Smith
1154ceaaa498SBarry Smith .seealso: `SNES, `SNESNCG`, `SNESNCGSetType()`
1155ceaaa498SBarry Smith M*/
11560a844d1aSPeter Brune typedef enum {
11570de8b71eSPeter Brune SNES_NCG_FR = 0,
11580de8b71eSPeter Brune SNES_NCG_PRP = 1,
11590de8b71eSPeter Brune SNES_NCG_HS = 2,
11600de8b71eSPeter Brune SNES_NCG_DY = 3,
11619371c9d4SSatish Balay SNES_NCG_CD = 4
11629371c9d4SSatish Balay } SNESNCGType;
11636a6fc655SJed Brown PETSC_EXTERN const char *const SNESNCGTypes[];
11640a844d1aSPeter Brune
1165014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESNCGSetType(SNES, SNESNCGType);
116613a62661SPeter Brune
1167ceaaa498SBarry Smith /*MC
1168ceaaa498SBarry Smith SNESQNScaleType - the scaling type used by `SNESQN`
1169ceaaa498SBarry Smith
1170ceaaa498SBarry Smith Values:
1171ceaaa498SBarry Smith + `SNES_QN_SCALE_NONE` - don't scale the problem
1172ceaaa498SBarry Smith . `SNES_QN_SCALE_SCALAR` - use Shanno scaling
1173ceaaa498SBarry Smith . `SNES_QN_SCALE_DIAGONAL` - scale with a diagonalized BFGS formula (see Gilbert and Lemarechal 1989), available
1174ceaaa498SBarry Smith - `SNES_QN_SCALE_JACOBIAN` - scale by solving a linear system coming from the Jacobian you provided with `SNESSetJacobian()`
1175ceaaa498SBarry Smith computed at the first iteration of `SNESQN` and at ever restart.
1176ceaaa498SBarry Smith
1177ceaaa498SBarry Smith Options Database Key:
1178ceaaa498SBarry Smith . -snes_qn_scale_type <diagonal,none,scalar,jacobian> - Scaling type
1179ceaaa498SBarry Smith
118095bd0b28SBarry Smith Level: intermediate
118195bd0b28SBarry Smith
1182ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNRestartType`
1183ceaaa498SBarry Smith M*/
11849371c9d4SSatish Balay typedef enum {
11859371c9d4SSatish Balay SNES_QN_SCALE_DEFAULT = 0,
11861efc8c45SPeter Brune SNES_QN_SCALE_NONE = 1,
118792f76d53SAlp Dener SNES_QN_SCALE_SCALAR = 2,
118892f76d53SAlp Dener SNES_QN_SCALE_DIAGONAL = 3,
11899371c9d4SSatish Balay SNES_QN_SCALE_JACOBIAN = 4
11909371c9d4SSatish Balay } SNESQNScaleType;
11916a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNScaleTypes[];
1192ceaaa498SBarry Smith
1193ceaaa498SBarry Smith /*MC
1194ceaaa498SBarry Smith SNESQNRestartType - the restart approached used by `SNESQN`
1195ceaaa498SBarry Smith
1196ceaaa498SBarry Smith Values:
1197ceaaa498SBarry Smith + `SNES_QN_RESTART_NONE` - never restart
1198ceaaa498SBarry Smith . `SNES_QN_RESTART_POWELL` - restart based upon descent criteria
1199ceaaa498SBarry Smith - `SNES_QN_RESTART_PERIODIC` - restart after a fixed number of iterations
1200ceaaa498SBarry Smith
1201ceaaa498SBarry Smith Options Database Keys:
1202ceaaa498SBarry Smith + -snes_qn_restart_type <powell,periodic,none> - set the restart type
1203ceaaa498SBarry Smith - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
1204ceaaa498SBarry Smith
120595bd0b28SBarry Smith Level: intermediate
120695bd0b28SBarry Smith
1207ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNType`, `SNESQNSetType()`, `SNESQNSetRestartType()`, `SNESQNScaleType`
1208ceaaa498SBarry Smith M*/
12099371c9d4SSatish Balay typedef enum {
12109371c9d4SSatish Balay SNES_QN_RESTART_DEFAULT = 0,
12111efc8c45SPeter Brune SNES_QN_RESTART_NONE = 1,
12121efc8c45SPeter Brune SNES_QN_RESTART_POWELL = 2,
12139371c9d4SSatish Balay SNES_QN_RESTART_PERIODIC = 3
12149371c9d4SSatish Balay } SNESQNRestartType;
12156a6fc655SJed Brown PETSC_EXTERN const char *const SNESQNRestartTypes[];
1216ceaaa498SBarry Smith
1217ceaaa498SBarry Smith /*MC
1218ceaaa498SBarry Smith SNESQNType - the type used by `SNESQN`
1219ceaaa498SBarry Smith
1220ceaaa498SBarry Smith Values:
1221ceaaa498SBarry Smith + `SNES_QN_LBFGS` - LBFGS variant
1222ceaaa498SBarry Smith . `SNES_QN_BROYDEN` - Broyden variant
1223ceaaa498SBarry Smith - `SNES_QN_BADBROYDEN` - Bad Broyden variant
1224ceaaa498SBarry Smith
1225ceaaa498SBarry Smith Options Database Key:
1226ceaaa498SBarry Smith . -snes_qn_type <lbfgs,broyden,badbroyden> - quasi-Newton type
1227ceaaa498SBarry Smith
122895bd0b28SBarry Smith Level: intermediate
122995bd0b28SBarry Smith
1230ceaaa498SBarry Smith .seealso: `SNES, `SNESQN`, `SNESQNSetScaleType()`, `SNESQNSetType()`, `SNESQNScaleType`, `SNESQNRestartType`, `SNESQNSetRestartType()`
1231ceaaa498SBarry Smith M*/
12329371c9d4SSatish Balay typedef enum {
12339371c9d4SSatish Balay SNES_QN_LBFGS = 0,
12341efc8c45SPeter Brune SNES_QN_BROYDEN = 1,
12351efc8c45SPeter Brune SNES_QN_BADBROYDEN = 2
12361efc8c45SPeter Brune } SNESQNType;
12371efc8c45SPeter Brune PETSC_EXTERN const char *const SNESQNTypes[];
12380c777b0cSPeter Brune
12391efc8c45SPeter Brune PETSC_EXTERN PetscErrorCode SNESQNSetType(SNES, SNESQNType);
1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetScaleType(SNES, SNESQNScaleType);
1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESQNSetRestartType(SNES, SNESQNRestartType);
12420c777b0cSPeter Brune
1243e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetType(SNES, PCASMType *);
1244e0331734SPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetType(SNES, PCASMType);
1245ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomains(SNES, PetscInt *, SNES *[], VecScatter *[], VecScatter *[], VecScatter *[]);
1246ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMSetSubdomains(SNES, PetscInt, SNES[], VecScatter[], VecScatter[], VecScatter[]);
1247610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetDamping(SNES, PetscReal);
1248610116beSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMGetDamping(SNES, PetscReal *);
1249ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode SNESNASMGetSubdomainVecs(SNES, PetscInt *, Vec *[], Vec *[], Vec *[], Vec *[]);
1250d728fb7dSPeter Brune PETSC_EXTERN PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES, PetscBool);
1251d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetSNES(SNES, PetscInt, SNES *);
1252d2dc0b00SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMGetNumber(SNES, PetscInt *);
1253f10b3e88SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESNASMSetWeight(SNES, Vec);
12540c777b0cSPeter Brune
1255ad05f2b7SBarry Smith /*E
1256ad05f2b7SBarry Smith SNESCompositeType - Determines how two or more preconditioners are composed with the `SNESType` of `SNESCOMPOSITE`
1257ad05f2b7SBarry Smith
1258ad05f2b7SBarry Smith Values:
1259ad05f2b7SBarry Smith + `SNES_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together
1260ad05f2b7SBarry Smith . `SNES_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
1261ad05f2b7SBarry Smith computed after the previous preconditioner application
1262ad05f2b7SBarry Smith - `SNES_COMPOSITE_ADDITIVEOPTIMAL` - uses a linear combination of the solutions obtained with each preconditioner that approximately minimize the function
1263ad05f2b7SBarry Smith value at the new iteration.
1264ad05f2b7SBarry Smith
1265ad05f2b7SBarry Smith Level: beginner
1266ad05f2b7SBarry Smith
1267ad05f2b7SBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `PCCompositeType`
1268ad05f2b7SBarry Smith E*/
12699371c9d4SSatish Balay typedef enum {
12709371c9d4SSatish Balay SNES_COMPOSITE_ADDITIVE,
12719371c9d4SSatish Balay SNES_COMPOSITE_MULTIPLICATIVE,
12729371c9d4SSatish Balay SNES_COMPOSITE_ADDITIVEOPTIMAL
12739371c9d4SSatish Balay } SNESCompositeType;
1274eed5f15bSPeter Brune PETSC_EXTERN const char *const SNESCompositeTypes[];
1275eed5f15bSPeter Brune
1276eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetType(SNES, SNESCompositeType);
1277eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeAddSNES(SNES, SNESType);
1278eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeGetSNES(SNES, PetscInt, SNES *);
1279a6b47ab3SPatrick Farrell PETSC_EXTERN PetscErrorCode SNESCompositeGetNumber(SNES, PetscInt *);
12808f626970SPeter Brune PETSC_EXTERN PetscErrorCode SNESCompositeSetDamping(SNES, PetscInt, PetscReal);
1281eed5f15bSPeter Brune
128268e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetDiscretisationInfo(SNES, PetscInt, DM *, PetscInt *, PetscInt *, const PetscInt **, const PetscInt *, PetscInt, const PetscInt *, PetscInt, const PetscInt *);
1283*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESPatchSetComputeOperator(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Mat, IS, PetscInt, const PetscInt *, const PetscInt *, PetscCtx), PetscCtx);
1284*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESPatchSetComputeFunction(SNES, PetscErrorCode (*func)(PC, PetscInt, Vec, Vec, IS, PetscInt, const PetscInt *, const PetscInt *, PetscCtx), PetscCtx);
1285*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESPatchSetConstructType(SNES, PCPatchConstructType, PetscErrorCode (*func)(PC, PetscInt *, IS **, IS *, PetscCtx), PetscCtx);
128668e6d53eSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESPatchSetCellNumbering(SNES, PetscSection);
1287561742edSMatthew G. Knepley
1288a055c8caSBarry Smith /*E
1289a055c8caSBarry Smith SNESFASType - Determines the type of nonlinear multigrid method that is run.
1290a055c8caSBarry Smith
1291a055c8caSBarry Smith Values:
129287497f52SBarry Smith + `SNES_FAS_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `SNESFASSetCycles()`
129387497f52SBarry Smith . `SNES_FAS_ADDITIVE` - additive FAS cycle
129487497f52SBarry Smith . `SNES_FAS_FULL` - full FAS cycle
129587497f52SBarry Smith - `SNES_FAS_KASKADE` - Kaskade FAS cycle
1296a055c8caSBarry Smith
129716a05f60SBarry Smith Level: beginner
129816a05f60SBarry Smith
12991cc06b55SBarry Smith .seealso: [](ch_snes), `SNESFAS`, `PCMGSetType()`, `PCMGType`
1300a055c8caSBarry Smith E*/
13019371c9d4SSatish Balay typedef enum {
13029371c9d4SSatish Balay SNES_FAS_MULTIPLICATIVE,
13039371c9d4SSatish Balay SNES_FAS_ADDITIVE,
13049371c9d4SSatish Balay SNES_FAS_FULL,
13059371c9d4SSatish Balay SNES_FAS_KASKADE
13069371c9d4SSatish Balay } SNESFASType;
1307a055c8caSBarry Smith PETSC_EXTERN const char *const SNESFASTypes[];
1308a055c8caSBarry Smith
1309a055c8caSBarry Smith /* called on the finest level FAS instance*/
1310a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetType(SNES, SNESFASType);
1311a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetType(SNES, SNESFASType *);
1312a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLevels(SNES, PetscInt, MPI_Comm *);
1313a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetLevels(SNES, PetscInt *);
1314a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCycleSNES(SNES, PetscInt, SNES *);
1315a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothUp(SNES, PetscInt);
1316a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetNumberSmoothDown(SNES, PetscInt);
1317a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetCycles(SNES, PetscInt);
1318d142ab34SLawrence Mitchell PETSC_EXTERN PetscErrorCode SNESFASSetMonitor(SNES, PetscViewerAndFormat *, PetscBool);
1319a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetLog(SNES, PetscBool);
1320a055c8caSBarry Smith
1321a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetGalerkin(SNES, PetscBool);
1322a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetGalerkin(SNES, PetscBool *);
1323*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGalerkinFunctionDefault(SNES, Vec, Vec, PetscCtx);
1324a055c8caSBarry Smith
1325a055c8caSBarry Smith /* called on any level -- "Cycle" FAS instance */
1326a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmoother(SNES, SNES *);
1327a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherUp(SNES, SNES *);
1328a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetSmootherDown(SNES, SNES *);
1329a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetCorrection(SNES, SNES *);
1330a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInterpolation(SNES, Mat *);
1331a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRestriction(SNES, Mat *);
1332a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetInjection(SNES, Mat *);
1333a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleGetRScale(SNES, Vec *);
1334a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleSetCycles(SNES, PetscInt);
1335a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCycleIsFine(SNES, PetscBool *);
1336a055c8caSBarry Smith
1337a055c8caSBarry Smith /* called on the (outer) finest level FAS to set/get parameters on any level instance */
1338a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInterpolation(SNES, PetscInt, Mat);
1339a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInterpolation(SNES, PetscInt, Mat *);
1340a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRestriction(SNES, PetscInt, Mat);
1341a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRestriction(SNES, PetscInt, Mat *);
1342a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetInjection(SNES, PetscInt, Mat);
1343a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetInjection(SNES, PetscInt, Mat *);
1344a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASSetRScale(SNES, PetscInt, Vec);
1345a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetRScale(SNES, PetscInt, Vec *);
1346020631bcSPeter Brune PETSC_EXTERN PetscErrorCode SNESFASSetContinuation(SNES, PetscBool);
1347a055c8caSBarry Smith
1348a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmoother(SNES, PetscInt, SNES *);
1349a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherUp(SNES, PetscInt, SNES *);
1350a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetSmootherDown(SNES, PetscInt, SNES *);
1351a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASGetCoarseSolve(SNES, SNES *);
1352a055c8caSBarry Smith
1353a055c8caSBarry Smith /* parameters for full FAS */
1354a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASFullSetDownSweep(SNES, PetscBool);
1355a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASCreateCoarseVec(SNES, Vec *);
1356a055c8caSBarry Smith PETSC_EXTERN PetscErrorCode SNESFASRestrict(SNES, Vec, Vec);
13576dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullSetTotal(SNES, PetscBool);
13586dfbafefSToby Isaac PETSC_EXTERN PetscErrorCode SNESFASFullGetTotal(SNES, PetscBool *);
1359a055c8caSBarry Smith
13602eace19aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetSNESVariableBounds(DM, SNES);
1361f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckDiscretization(SNES, DM, PetscReal, Vec, PetscReal, PetscReal[]);
1362f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckResidual(SNES, DM, Vec, PetscReal, PetscReal *);
1363f3c94560SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckJacobian(SNES, DM, Vec, PetscReal, PetscBool *, PetscReal *);
13647f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSNESCheckFromOptions(SNES, Vec);
1365*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESComputeJacobianAction(DM, Vec, Vec, Vec, PetscCtx);
1366*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSNESCreateJacobianMF(DM, Vec, PetscCtx, Mat *);
136797276fddSZach Atkins
1368*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonALSetFunction(SNES, SNESFunctionFn *, PetscCtx);
1369*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode SNESNewtonALGetFunction(SNES, SNESFunctionFn **, PetscCtxRt);
137097276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALComputeFunction(SNES, Vec, Vec);
137197276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALGetLoadParameter(SNES, PetscReal *);
137297276fddSZach Atkins
137397276fddSZach Atkins /*MC
137497276fddSZach Atkins SNESNewtonALCorrectionType - the approach used by `SNESNEWTONAL` to determine
137597276fddSZach Atkins the correction to the current increment. While the exact correction satisfies
137697276fddSZach Atkins the constraint surface at every iteration, it also requires solving a quadratic
137797276fddSZach Atkins equation which may not have real roots. Conversely, the normal correction is more
137897276fddSZach Atkins efficient and always yields a real correction and is the default.
137997276fddSZach Atkins
138097276fddSZach Atkins Values:
138197276fddSZach Atkins + `SNES_NEWTONAL_CORRECTION_EXACT` - choose the correction which exactly satisfies the constraint
1382d7c1f440SPierre Jolivet - `SNES_NEWTONAL_CORRECTION_NORMAL` - choose the correction in the updated normal hyper-surface to the constraint surface
138397276fddSZach Atkins
138497276fddSZach Atkins Options Database Key:
138597276fddSZach Atkins . -snes_newtonal_correction_type <exact> - select type from <exact,normal>
138697276fddSZach Atkins
138797276fddSZach Atkins Level: intermediate
138897276fddSZach Atkins
138997276fddSZach Atkins .seealso: `SNES`, `SNESNEWTONAL`, `SNESNewtonALSetCorrectionType()`
139097276fddSZach Atkins M*/
139197276fddSZach Atkins typedef enum {
139297276fddSZach Atkins SNES_NEWTONAL_CORRECTION_EXACT = 0,
139397276fddSZach Atkins SNES_NEWTONAL_CORRECTION_NORMAL = 1,
139497276fddSZach Atkins } SNESNewtonALCorrectionType;
139597276fddSZach Atkins PETSC_EXTERN const char *const SNESNewtonALCorrectionTypes[];
139697276fddSZach Atkins
139797276fddSZach Atkins PETSC_EXTERN PetscErrorCode SNESNewtonALSetCorrectionType(SNES, SNESNewtonALCorrectionType);
1398