xref: /petsc/src/snes/impls/gs/snesgs.c (revision bd89dbf26d8a5efecb980364933175da61864cd7)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h"  I*/
26cd56b8bSPeter Brune 
3b6266c6eSPeter Brune /*@
4f6dfbefdSBarry Smith   SNESNGSSetTolerances - Sets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG`
5b6266c6eSPeter Brune 
6c3339decSBarry Smith   Logically Collective
7b6266c6eSPeter Brune 
8b6266c6eSPeter Brune   Input Parameters:
9f6dfbefdSBarry Smith + snes   - the `SNES` context
10b6266c6eSPeter Brune . abstol - absolute convergence tolerance
11b6266c6eSPeter Brune . rtol   - relative convergence tolerance
12b6266c6eSPeter Brune . stol   - convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
13b6266c6eSPeter Brune - maxit  - maximum number of iterations
14b6266c6eSPeter Brune 
15b6266c6eSPeter Brune   Options Database Keys:
16be95d8f1SBarry Smith + -snes_ngs_atol <abstol> - Sets abstol
17be95d8f1SBarry Smith . -snes_ngs_rtol <rtol>   - Sets rtol
18be95d8f1SBarry Smith . -snes_ngs_stol <stol>   - Sets stol
19b6266c6eSPeter Brune - -snes_max_it <maxit>    - Sets maxit
20b6266c6eSPeter Brune 
21b6266c6eSPeter Brune   Level: intermediate
22b6266c6eSPeter Brune 
2377e5a1f9SBarry Smith   Notes:
2477e5a1f9SBarry Smith   Use `PETSC_CURRENT` to retain the value for any parameter
2577e5a1f9SBarry Smith 
2677e5a1f9SBarry Smith   All parameters must be non-negative
2777e5a1f9SBarry Smith 
2877e5a1f9SBarry Smith   Developer Note:
2977e5a1f9SBarry Smith   Why can't the values set with `SNESSetTolerances()` be used?
3077e5a1f9SBarry Smith 
313201ab8dSStefano Zampini .seealso: [](ch_snes), `SNES`, `SNESNCG`
32b6266c6eSPeter Brune @*/
SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)33d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSSetTolerances(SNES snes, PetscReal abstol, PetscReal rtol, PetscReal stol, PetscInt maxit)
34d71ae5a4SJacob Faibussowitsch {
35be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
36b6266c6eSPeter Brune 
37b6266c6eSPeter Brune   PetscFunctionBegin;
38b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
39b6266c6eSPeter Brune 
4077e5a1f9SBarry Smith   if (abstol != (PetscReal)PETSC_CURRENT) {
4108401ef6SPierre Jolivet     PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
42b6266c6eSPeter Brune     gs->abstol = abstol;
43b6266c6eSPeter Brune   }
4477e5a1f9SBarry Smith   if (rtol != (PetscReal)PETSC_CURRENT) {
450b121fc5SBarry Smith     PetscCheck(rtol >= 0.0 && rtol < 1.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %g must be non-negative and less than 1.0", (double)rtol);
46b6266c6eSPeter Brune     gs->rtol = rtol;
47b6266c6eSPeter Brune   }
4877e5a1f9SBarry Smith   if (stol != (PetscReal)PETSC_CURRENT) {
4908401ef6SPierre Jolivet     PetscCheck(stol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Step tolerance %g must be non-negative", (double)stol);
50b6266c6eSPeter Brune     gs->stol = stol;
51b6266c6eSPeter Brune   }
5277e5a1f9SBarry Smith   if (maxit != PETSC_CURRENT) {
5363a3b9bcSJacob Faibussowitsch     PetscCheck(maxit >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxit);
54b6266c6eSPeter Brune     gs->max_its = maxit;
55b6266c6eSPeter Brune   }
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57b6266c6eSPeter Brune }
58b6266c6eSPeter Brune 
59b6266c6eSPeter Brune /*@
60f6dfbefdSBarry Smith   SNESNGSGetTolerances - Gets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG`
61b6266c6eSPeter Brune 
62b6266c6eSPeter Brune   Not Collective
63b6266c6eSPeter Brune 
64b6266c6eSPeter Brune   Input Parameters:
65f6dfbefdSBarry Smith + snes  - the `SNES` context
66b6266c6eSPeter Brune . atol  - absolute convergence tolerance
67b6266c6eSPeter Brune . rtol  - relative convergence tolerance
68b6266c6eSPeter Brune . stol  - convergence tolerance in terms of the norm
69b6266c6eSPeter Brune            of the change in the solution between steps
70b6266c6eSPeter Brune - maxit - maximum number of iterations
71b6266c6eSPeter Brune 
722fe279fdSBarry Smith   Level: intermediate
732fe279fdSBarry Smith 
74f6dfbefdSBarry Smith   Note:
75420bcc1bSBarry Smith   The user can specify `NULL` for any parameter that is not needed.
76b6266c6eSPeter Brune 
77420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetTolerances()`
78b6266c6eSPeter Brune @*/
SNESNGSGetTolerances(SNES snes,PetscReal * atol,PetscReal * rtol,PetscReal * stol,PetscInt * maxit)79d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSGetTolerances(SNES snes, PetscReal *atol, PetscReal *rtol, PetscReal *stol, PetscInt *maxit)
80d71ae5a4SJacob Faibussowitsch {
81be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
82b6266c6eSPeter Brune 
83b6266c6eSPeter Brune   PetscFunctionBegin;
84b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
85b6266c6eSPeter Brune   if (atol) *atol = gs->abstol;
86b6266c6eSPeter Brune   if (rtol) *rtol = gs->rtol;
87b6266c6eSPeter Brune   if (stol) *stol = gs->stol;
88b6266c6eSPeter Brune   if (maxit) *maxit = gs->max_its;
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90b6266c6eSPeter Brune }
916cd56b8bSPeter Brune 
926cd56b8bSPeter Brune /*@
93f6dfbefdSBarry Smith   SNESNGSSetSweeps - Sets the number of sweeps of nonlinear GS to use in `SNESNCG`
946cd56b8bSPeter Brune 
95420bcc1bSBarry Smith   Logically Collective
96420bcc1bSBarry Smith 
976cd56b8bSPeter Brune   Input Parameters:
98f6dfbefdSBarry Smith + snes   - the `SNES` context
99f6dfbefdSBarry Smith - sweeps - the number of sweeps of nonlinear GS to perform.
100f6dfbefdSBarry Smith 
101f6dfbefdSBarry Smith   Options Database Key:
102f6dfbefdSBarry Smith . -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply
1036cd56b8bSPeter Brune 
1046cd56b8bSPeter Brune   Level: intermediate
1056cd56b8bSPeter Brune 
106420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()`
1076cd56b8bSPeter Brune @*/
SNESNGSSetSweeps(SNES snes,PetscInt sweeps)108d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
109d71ae5a4SJacob Faibussowitsch {
110be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1116cd56b8bSPeter Brune 
1126cd56b8bSPeter Brune   PetscFunctionBegin;
1136cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1146cd56b8bSPeter Brune   gs->sweeps = sweeps;
1153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1166cd56b8bSPeter Brune }
1176cd56b8bSPeter Brune 
1186cd56b8bSPeter Brune /*@
119f6dfbefdSBarry Smith   SNESNGSGetSweeps - Gets the number of sweeps nonlinear GS will use in `SNESNCG`
1206cd56b8bSPeter Brune 
1212fe279fdSBarry Smith   Input Parameter:
122f6dfbefdSBarry Smith . snes - the `SNES` context
1236cd56b8bSPeter Brune 
1242fe279fdSBarry Smith   Output Parameter:
125f6dfbefdSBarry Smith . sweeps - the number of sweeps of nonlinear GS to perform.
1266cd56b8bSPeter Brune 
1276cd56b8bSPeter Brune   Level: intermediate
1286cd56b8bSPeter Brune 
129420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()`
1306cd56b8bSPeter Brune @*/
SNESNGSGetSweeps(SNES snes,PetscInt * sweeps)131d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps)
132d71ae5a4SJacob Faibussowitsch {
133be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1346cd56b8bSPeter Brune 
1356cd56b8bSPeter Brune   PetscFunctionBegin;
1366cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1376cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1396cd56b8bSPeter Brune }
1406cd56b8bSPeter Brune 
SNESReset_NGS(SNES snes)14166976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NGS(SNES snes)
142d71ae5a4SJacob Faibussowitsch {
1438a86d3c5SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1448a86d3c5SBarry Smith 
1453542a6bcSPeter Brune   PetscFunctionBegin;
1469566063dSJacob Faibussowitsch   PetscCall(ISColoringDestroy(&gs->coloring));
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1483542a6bcSPeter Brune }
1493542a6bcSPeter Brune 
SNESDestroy_NGS(SNES snes)15066976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NGS(SNES snes)
151d71ae5a4SJacob Faibussowitsch {
1523542a6bcSPeter Brune   PetscFunctionBegin;
1539566063dSJacob Faibussowitsch   PetscCall(SNESReset_NGS(snes));
1549566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
1553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1563542a6bcSPeter Brune }
1573542a6bcSPeter Brune 
SNESSetUp_NGS(SNES snes)15866976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NGS(SNES snes)
159d71ae5a4SJacob Faibussowitsch {
160ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
161ef107cf5SPeter Brune 
1623542a6bcSPeter Brune   PetscFunctionBegin;
1639566063dSJacob Faibussowitsch   PetscCall(SNESGetNGS(snes, &f, NULL));
16448a46eb9SPierre Jolivet   if (!f) PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1663542a6bcSPeter Brune }
1673542a6bcSPeter Brune 
SNESSetFromOptions_NGS(SNES snes,PetscOptionItems PetscOptionsObject)168ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems PetscOptionsObject)
169d71ae5a4SJacob Faibussowitsch {
170be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
17177e5a1f9SBarry Smith   PetscInt  sweeps, max_its = PETSC_CURRENT;
17277e5a1f9SBarry Smith   PetscReal rtol = PETSC_CURRENT, atol = PETSC_CURRENT, stol = PETSC_CURRENT;
173b6266c6eSPeter Brune   PetscBool flg, flg1, flg2, flg3;
1743542a6bcSPeter Brune 
1753542a6bcSPeter Brune   PetscFunctionBegin;
176d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options");
1776cd56b8bSPeter Brune   /* GS Options */
1789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg));
1791baa6e33SBarry Smith   if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps));
1809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg));
1819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1));
1829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2));
1839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3));
18448a46eb9SPierre Jolivet   if (flg || flg1 || flg2 || flg3) PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its));
185ef107cf5SPeter Brune   flg = PETSC_FALSE;
1869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL));
187ef107cf5SPeter Brune   if (flg) {
1889566063dSJacob Faibussowitsch     PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1899566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n"));
190ef107cf5SPeter Brune   }
1919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL));
1929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngs_secant_mat_coloring", "Use the graph coloring of the Jacobian for the secant GS", "", gs->secant_mat, &gs->secant_mat, &flg));
193ef107cf5SPeter Brune 
194d0609cedSBarry Smith   PetscOptionsHeadEnd();
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1963542a6bcSPeter Brune }
1973542a6bcSPeter Brune 
SNESView_NGS(SNES snes,PetscViewer viewer)19866976f2fSJacob Faibussowitsch static PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
199d71ae5a4SJacob Faibussowitsch {
20098272f38SBarry Smith   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
20198272f38SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
2029f196a02SMartin Diehl   PetscBool isascii;
20398272f38SBarry Smith 
2043542a6bcSPeter Brune   PetscFunctionBegin;
2059f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2069f196a02SMartin Diehl   if (isascii) {
2079566063dSJacob Faibussowitsch     PetscCall(DMSNESGetNGS(snes->dm, &f, NULL));
20848a46eb9SPierre Jolivet     if (f == SNESComputeNGSDefaultSecant) PetscCall(PetscViewerASCIIPrintf(viewer, "  Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h));
20998272f38SBarry Smith   }
2103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2113542a6bcSPeter Brune }
2123542a6bcSPeter Brune 
SNESSolve_NGS(SNES snes)21366976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NGS(SNES snes)
214d71ae5a4SJacob Faibussowitsch {
2153542a6bcSPeter Brune   Vec              F;
2163542a6bcSPeter Brune   Vec              X;
2173542a6bcSPeter Brune   Vec              B;
2183542a6bcSPeter Brune   PetscInt         i;
21906e07b1aSPeter Brune   PetscReal        fnorm;
220365a6726SPeter Brune   SNESNormSchedule normschedule;
2213542a6bcSPeter Brune 
2223542a6bcSPeter Brune   PetscFunctionBegin;
2230b121fc5SBarry Smith   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
224c579b300SPatrick Farrell 
2259566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
2263542a6bcSPeter Brune   X = snes->vec_sol;
2273542a6bcSPeter Brune   F = snes->vec_func;
2283542a6bcSPeter Brune   B = snes->vec_rhs;
2293542a6bcSPeter Brune 
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
231534ebe21SPeter Brune   snes->iter = 0;
232534ebe21SPeter Brune   snes->norm = 0.;
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
234526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
235534ebe21SPeter Brune 
2369566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
237365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2383542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
239ac530a7eSPierre Jolivet     if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
240ac530a7eSPierre Jolivet     else snes->vec_func_init_set = PETSC_FALSE;
241e4ed7901SPeter Brune 
2429566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
243*76c63389SBarry Smith     SNESCheckFunctionDomainError(snes, fnorm);
2449566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
245e4ed7901SPeter Brune     snes->iter = 0;
2463542a6bcSPeter Brune     snes->norm = fnorm;
2479566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2489566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
2493542a6bcSPeter Brune 
2503542a6bcSPeter Brune     /* test convergence */
251d76a863bSStefano Zampini     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
2522d157150SStefano Zampini     PetscCall(SNESMonitor(snes, 0, snes->norm));
2533ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
254534ebe21SPeter Brune   } else {
2559566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2569566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
257534ebe21SPeter Brune   }
258534ebe21SPeter Brune 
2593542a6bcSPeter Brune   /* Call general purpose update function */
260dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
261534ebe21SPeter Brune 
2623542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
2639566063dSJacob Faibussowitsch     PetscCall(SNESComputeNGS(snes, B, X));
2644b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
265365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2669566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2679566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
268*76c63389SBarry Smith       SNESCheckFunctionDomainError(snes, fnorm);
2692d157150SStefano Zampini     }
2703542a6bcSPeter Brune     /* Monitor convergence */
2719566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
2723542a6bcSPeter Brune     snes->iter = i + 1;
2733542a6bcSPeter Brune     snes->norm = fnorm;
2749566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2752d157150SStefano Zampini     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
2763542a6bcSPeter Brune     /* Test for convergence */
2772d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
2782d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2793ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
2803542a6bcSPeter Brune     /* Call general purpose update function */
281dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
2823542a6bcSPeter Brune   }
2833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2843542a6bcSPeter Brune }
2853542a6bcSPeter Brune 
2863542a6bcSPeter Brune /*MC
287420bcc1bSBarry Smith   SNESNGS - Either calls the user-provided Gauss-Seidel solution routine provided with `SNESSetNGS()` or does a finite difference secant approximation
2881d27aa22SBarry Smith             using coloring {cite}`bruneknepleysmithtu15`.
2893542a6bcSPeter Brune 
2903542a6bcSPeter Brune    Level: advanced
2913542a6bcSPeter Brune 
292f6dfbefdSBarry Smith   Options Database Keys:
293f6dfbefdSBarry Smith +   -snes_ngs_sweeps <n>                                                          - Number of sweeps of nonlinear GS to apply
294f6dfbefdSBarry Smith .    -snes_ngs_atol <atol>                                                        - Absolute residual tolerance for nonlinear GS iteration
295f6dfbefdSBarry Smith .    -snes_ngs_rtol <rtol>                                                        - Relative residual tolerance for nonlinear GS iteration
296f6dfbefdSBarry Smith .    -snes_ngs_stol <stol>                                                        - Absolute update tolerance for nonlinear GS iteration
297f6dfbefdSBarry Smith .    -snes_ngs_max_it <maxit>                                                     - Maximum number of sweeps of nonlinea GS to apply
2981d27aa22SBarry Smith .    -snes_ngs_secant                                                             - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine,
2991d27aa22SBarry Smith                                                                                     this is used by default if no user provided Gauss-Seidel routine is available.
3001d27aa22SBarry Smith                                                                                     Requires either that a `DM` that can compute a coloring
30198272f38SBarry Smith                                                                                     is available or a Jacobian sparse matrix is provided (from which to get the coloring).
30278e24b28SBarry Smith .    -snes_ngs_secant_h <h>                                                       - Differencing parameter for secant approximation
3031d27aa22SBarry Smith .    -snes_ngs_secant_mat_coloring                                                - Use the graph coloring of the Jacobian for the secant GS even if a `DM` is available.
304a5b23f4aSJose E. Roman -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed
30598272f38SBarry Smith 
3063542a6bcSPeter Brune   Notes:
307f6dfbefdSBarry Smith   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with `SNESGetNPC()`, it will have
3083542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
3093542a6bcSPeter Brune 
310f6dfbefdSBarry Smith   By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with `SNESSetNormSchedule()`
31177e5a1f9SBarry Smith   or `-snes_norm_schedule none`
312f6dfbefdSBarry Smith 
313420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`,
314420bcc1bSBarry Smith           `SNESSetNormSchedule()`, `SNESNGSGetTolerances()`, `SNESNGSSetSweeps()`
3153542a6bcSPeter Brune M*/
3163542a6bcSPeter Brune 
SNESCreate_NGS(SNES snes)317d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
318d71ae5a4SJacob Faibussowitsch {
319be95d8f1SBarry Smith   SNES_NGS *gs;
3203542a6bcSPeter Brune 
3213542a6bcSPeter Brune   PetscFunctionBegin;
322be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
323be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
324be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
325be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
326be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
327be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3283542a6bcSPeter Brune 
3293542a6bcSPeter Brune   snes->usesksp                     = PETSC_FALSE;
330efd4aadfSBarry Smith   snes->usesnpc                     = PETSC_FALSE;
3314fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
3324fc747eaSLawrence Mitchell 
33377e5a1f9SBarry Smith   PetscCall(SNESParametersInitialize(snes));
33477e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_funcs, 10000);
33577e5a1f9SBarry Smith   PetscObjectParameterSetDefault(snes, max_its, 10000);
3360e444f03SPeter Brune 
3374dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&gs));
3386cd56b8bSPeter Brune   gs->sweeps  = 1;
339b6266c6eSPeter Brune   gs->rtol    = 1e-5;
34078e24b28SBarry Smith   gs->abstol  = PETSC_MACHINE_EPSILON;
34178e24b28SBarry Smith   gs->stol    = 1000 * PETSC_MACHINE_EPSILON;
342b6266c6eSPeter Brune   gs->max_its = 50;
34378e24b28SBarry Smith   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
3443542a6bcSPeter Brune   snes->data  = (void *)gs;
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3463542a6bcSPeter Brune }
347