xref: /petsc/src/snes/impls/gs/snesgs.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
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 
23420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetTrustRegionTolerance()`
24b6266c6eSPeter Brune @*/
25d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSSetTolerances(SNES snes, PetscReal abstol, PetscReal rtol, PetscReal stol, PetscInt maxit)
26d71ae5a4SJacob Faibussowitsch {
27be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
28b6266c6eSPeter Brune 
29b6266c6eSPeter Brune   PetscFunctionBegin;
30b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
31b6266c6eSPeter Brune 
3213bcc0bdSJacob Faibussowitsch   if (abstol != (PetscReal)PETSC_DEFAULT) {
3308401ef6SPierre Jolivet     PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
34b6266c6eSPeter Brune     gs->abstol = abstol;
35b6266c6eSPeter Brune   }
3613bcc0bdSJacob Faibussowitsch   if (rtol != (PetscReal)PETSC_DEFAULT) {
370b121fc5SBarry 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);
38b6266c6eSPeter Brune     gs->rtol = rtol;
39b6266c6eSPeter Brune   }
4013bcc0bdSJacob Faibussowitsch   if (stol != (PetscReal)PETSC_DEFAULT) {
4108401ef6SPierre Jolivet     PetscCheck(stol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Step tolerance %g must be non-negative", (double)stol);
42b6266c6eSPeter Brune     gs->stol = stol;
43b6266c6eSPeter Brune   }
44b6266c6eSPeter Brune   if (maxit != PETSC_DEFAULT) {
4563a3b9bcSJacob Faibussowitsch     PetscCheck(maxit >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxit);
46b6266c6eSPeter Brune     gs->max_its = maxit;
47b6266c6eSPeter Brune   }
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49b6266c6eSPeter Brune }
50b6266c6eSPeter Brune 
51b6266c6eSPeter Brune /*@
52f6dfbefdSBarry Smith   SNESNGSGetTolerances - Gets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG`
53b6266c6eSPeter Brune 
54b6266c6eSPeter Brune   Not Collective
55b6266c6eSPeter Brune 
56b6266c6eSPeter Brune   Input Parameters:
57f6dfbefdSBarry Smith + snes  - the `SNES` context
58b6266c6eSPeter Brune . atol  - absolute convergence tolerance
59b6266c6eSPeter Brune . rtol  - relative convergence tolerance
60b6266c6eSPeter Brune . stol  - convergence tolerance in terms of the norm
61b6266c6eSPeter Brune            of the change in the solution between steps
62b6266c6eSPeter Brune - maxit - maximum number of iterations
63b6266c6eSPeter Brune 
642fe279fdSBarry Smith   Level: intermediate
652fe279fdSBarry Smith 
66f6dfbefdSBarry Smith   Note:
67420bcc1bSBarry Smith   The user can specify `NULL` for any parameter that is not needed.
68b6266c6eSPeter Brune 
69420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetTolerances()`
70b6266c6eSPeter Brune @*/
71d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSGetTolerances(SNES snes, PetscReal *atol, PetscReal *rtol, PetscReal *stol, PetscInt *maxit)
72d71ae5a4SJacob Faibussowitsch {
73be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
74b6266c6eSPeter Brune 
75b6266c6eSPeter Brune   PetscFunctionBegin;
76b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
77b6266c6eSPeter Brune   if (atol) *atol = gs->abstol;
78b6266c6eSPeter Brune   if (rtol) *rtol = gs->rtol;
79b6266c6eSPeter Brune   if (stol) *stol = gs->stol;
80b6266c6eSPeter Brune   if (maxit) *maxit = gs->max_its;
813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
82b6266c6eSPeter Brune }
836cd56b8bSPeter Brune 
846cd56b8bSPeter Brune /*@
85f6dfbefdSBarry Smith   SNESNGSSetSweeps - Sets the number of sweeps of nonlinear GS to use in `SNESNCG`
866cd56b8bSPeter Brune 
87420bcc1bSBarry Smith   Logically Collective
88420bcc1bSBarry Smith 
896cd56b8bSPeter Brune   Input Parameters:
90f6dfbefdSBarry Smith + snes   - the `SNES` context
91f6dfbefdSBarry Smith - sweeps - the number of sweeps of nonlinear GS to perform.
92f6dfbefdSBarry Smith 
93f6dfbefdSBarry Smith   Options Database Key:
94f6dfbefdSBarry Smith . -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply
956cd56b8bSPeter Brune 
966cd56b8bSPeter Brune   Level: intermediate
976cd56b8bSPeter Brune 
98420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()`
996cd56b8bSPeter Brune @*/
100d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
101d71ae5a4SJacob Faibussowitsch {
102be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1036cd56b8bSPeter Brune 
1046cd56b8bSPeter Brune   PetscFunctionBegin;
1056cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1066cd56b8bSPeter Brune   gs->sweeps = sweeps;
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1086cd56b8bSPeter Brune }
1096cd56b8bSPeter Brune 
1106cd56b8bSPeter Brune /*@
111f6dfbefdSBarry Smith   SNESNGSGetSweeps - Gets the number of sweeps nonlinear GS will use in `SNESNCG`
1126cd56b8bSPeter Brune 
1132fe279fdSBarry Smith   Input Parameter:
114f6dfbefdSBarry Smith . snes - the `SNES` context
1156cd56b8bSPeter Brune 
1162fe279fdSBarry Smith   Output Parameter:
117f6dfbefdSBarry Smith . sweeps - the number of sweeps of nonlinear GS to perform.
1186cd56b8bSPeter Brune 
1196cd56b8bSPeter Brune   Level: intermediate
1206cd56b8bSPeter Brune 
121420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()`
1226cd56b8bSPeter Brune @*/
123d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps)
124d71ae5a4SJacob Faibussowitsch {
125be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1266cd56b8bSPeter Brune 
1276cd56b8bSPeter Brune   PetscFunctionBegin;
1286cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1296cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1316cd56b8bSPeter Brune }
1326cd56b8bSPeter Brune 
13366976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NGS(SNES snes)
134d71ae5a4SJacob Faibussowitsch {
1358a86d3c5SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1368a86d3c5SBarry Smith 
1373542a6bcSPeter Brune   PetscFunctionBegin;
1389566063dSJacob Faibussowitsch   PetscCall(ISColoringDestroy(&gs->coloring));
1393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1403542a6bcSPeter Brune }
1413542a6bcSPeter Brune 
14266976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NGS(SNES snes)
143d71ae5a4SJacob Faibussowitsch {
1443542a6bcSPeter Brune   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(SNESReset_NGS(snes));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1483542a6bcSPeter Brune }
1493542a6bcSPeter Brune 
15066976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NGS(SNES snes)
151d71ae5a4SJacob Faibussowitsch {
152ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
153ef107cf5SPeter Brune 
1543542a6bcSPeter Brune   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(SNESGetNGS(snes, &f, NULL));
15648a46eb9SPierre Jolivet   if (!f) PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1583542a6bcSPeter Brune }
1593542a6bcSPeter Brune 
16066976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems *PetscOptionsObject)
161d71ae5a4SJacob Faibussowitsch {
162be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
163b6266c6eSPeter Brune   PetscInt  sweeps, max_its = PETSC_DEFAULT;
164b6266c6eSPeter Brune   PetscReal rtol = PETSC_DEFAULT, atol = PETSC_DEFAULT, stol = PETSC_DEFAULT;
165b6266c6eSPeter Brune   PetscBool flg, flg1, flg2, flg3;
1663542a6bcSPeter Brune 
1673542a6bcSPeter Brune   PetscFunctionBegin;
168d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options");
1696cd56b8bSPeter Brune   /* GS Options */
1709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg));
1711baa6e33SBarry Smith   if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps));
1729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg));
1739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1));
1749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2));
1759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3));
17648a46eb9SPierre Jolivet   if (flg || flg1 || flg2 || flg3) PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its));
177ef107cf5SPeter Brune   flg = PETSC_FALSE;
1789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL));
179ef107cf5SPeter Brune   if (flg) {
1809566063dSJacob Faibussowitsch     PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1819566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n"));
182ef107cf5SPeter Brune   }
1839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL));
1849566063dSJacob 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));
185ef107cf5SPeter Brune 
186d0609cedSBarry Smith   PetscOptionsHeadEnd();
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1883542a6bcSPeter Brune }
1893542a6bcSPeter Brune 
19066976f2fSJacob Faibussowitsch static PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
191d71ae5a4SJacob Faibussowitsch {
19298272f38SBarry Smith   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
19398272f38SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
19498272f38SBarry Smith   PetscBool iascii;
19598272f38SBarry Smith 
1963542a6bcSPeter Brune   PetscFunctionBegin;
1979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
19898272f38SBarry Smith   if (iascii) {
1999566063dSJacob Faibussowitsch     PetscCall(DMSNESGetNGS(snes->dm, &f, NULL));
20048a46eb9SPierre Jolivet     if (f == SNESComputeNGSDefaultSecant) PetscCall(PetscViewerASCIIPrintf(viewer, "  Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h));
20198272f38SBarry Smith   }
2023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2033542a6bcSPeter Brune }
2043542a6bcSPeter Brune 
20566976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NGS(SNES snes)
206d71ae5a4SJacob Faibussowitsch {
2073542a6bcSPeter Brune   Vec              F;
2083542a6bcSPeter Brune   Vec              X;
2093542a6bcSPeter Brune   Vec              B;
2103542a6bcSPeter Brune   PetscInt         i;
21106e07b1aSPeter Brune   PetscReal        fnorm;
212365a6726SPeter Brune   SNESNormSchedule normschedule;
2133542a6bcSPeter Brune 
2143542a6bcSPeter Brune   PetscFunctionBegin;
215c579b300SPatrick Farrell 
2160b121fc5SBarry 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);
217c579b300SPatrick Farrell 
2189566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
2193542a6bcSPeter Brune   X = snes->vec_sol;
2203542a6bcSPeter Brune   F = snes->vec_func;
2213542a6bcSPeter Brune   B = snes->vec_rhs;
2223542a6bcSPeter Brune 
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
224534ebe21SPeter Brune   snes->iter = 0;
225534ebe21SPeter Brune   snes->norm = 0.;
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
227526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
228534ebe21SPeter Brune 
2299566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
230365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2313542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
232e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2339566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2341aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
235e4ed7901SPeter Brune 
2369566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
237422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
239e4ed7901SPeter Brune     snes->iter = 0;
2403542a6bcSPeter Brune     snes->norm = fnorm;
2419566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2429566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
2433542a6bcSPeter Brune 
2443542a6bcSPeter Brune     /* test convergence */
245d76a863bSStefano Zampini     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
2462d157150SStefano Zampini     PetscCall(SNESMonitor(snes, 0, snes->norm));
2473ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
248534ebe21SPeter Brune   } else {
2499566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2509566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
251534ebe21SPeter Brune   }
252534ebe21SPeter Brune 
2533542a6bcSPeter Brune   /* Call general purpose update function */
254dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
255534ebe21SPeter Brune 
2563542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
2579566063dSJacob Faibussowitsch     PetscCall(SNESComputeNGS(snes, B, X));
2584b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
259365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2609566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2619566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
262422a814eSBarry Smith       SNESCheckFunctionNorm(snes, fnorm);
2632d157150SStefano Zampini     }
2643542a6bcSPeter Brune     /* Monitor convergence */
2659566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
2663542a6bcSPeter Brune     snes->iter = i + 1;
2673542a6bcSPeter Brune     snes->norm = fnorm;
2689566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2692d157150SStefano Zampini     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
2703542a6bcSPeter Brune     /* Test for convergence */
2712d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
2722d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2733ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
2743542a6bcSPeter Brune     /* Call general purpose update function */
275dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
2763542a6bcSPeter Brune   }
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2783542a6bcSPeter Brune }
2793542a6bcSPeter Brune 
2803542a6bcSPeter Brune /*MC
281420bcc1bSBarry Smith   SNESNGS - Either calls the user-provided Gauss-Seidel solution routine provided with `SNESSetNGS()` or does a finite difference secant approximation
282*1d27aa22SBarry Smith             using coloring {cite}`bruneknepleysmithtu15`.
2833542a6bcSPeter Brune 
2843542a6bcSPeter Brune    Level: advanced
2853542a6bcSPeter Brune 
286f6dfbefdSBarry Smith   Options Database Keys:
287f6dfbefdSBarry Smith +   -snes_ngs_sweeps <n>                                                          - Number of sweeps of nonlinear GS to apply
288f6dfbefdSBarry Smith .    -snes_ngs_atol <atol>                                                        - Absolute residual tolerance for nonlinear GS iteration
289f6dfbefdSBarry Smith .    -snes_ngs_rtol <rtol>                                                        - Relative residual tolerance for nonlinear GS iteration
290f6dfbefdSBarry Smith .    -snes_ngs_stol <stol>                                                        - Absolute update tolerance for nonlinear GS iteration
291f6dfbefdSBarry Smith .    -snes_ngs_max_it <maxit>                                                     - Maximum number of sweeps of nonlinea GS to apply
292*1d27aa22SBarry Smith .    -snes_ngs_secant                                                             - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine,
293*1d27aa22SBarry Smith                                                                                     this is used by default if no user provided Gauss-Seidel routine is available.
294*1d27aa22SBarry Smith                                                                                     Requires either that a `DM` that can compute a coloring
29598272f38SBarry Smith                                                                                     is available or a Jacobian sparse matrix is provided (from which to get the coloring).
29678e24b28SBarry Smith .    -snes_ngs_secant_h <h>                                                       - Differencing parameter for secant approximation
297*1d27aa22SBarry Smith .    -snes_ngs_secant_mat_coloring                                                - Use the graph coloring of the Jacobian for the secant GS even if a `DM` is available.
298a5b23f4aSJose E. Roman -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed
29998272f38SBarry Smith 
3003542a6bcSPeter Brune   Notes:
301f6dfbefdSBarry Smith   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with `SNESGetNPC()`, it will have
3023542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
3033542a6bcSPeter Brune 
304f6dfbefdSBarry Smith   By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with `SNESSetNormSchedule()`
305f6dfbefdSBarry Smith   or -snes_norm_schedule none
306f6dfbefdSBarry Smith 
307420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`,
308420bcc1bSBarry Smith           `SNESSetNormSchedule()`, `SNESNGSGetTolerances()`, `SNESNGSSetSweeps()`
3093542a6bcSPeter Brune M*/
3103542a6bcSPeter Brune 
311d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
312d71ae5a4SJacob Faibussowitsch {
313be95d8f1SBarry Smith   SNES_NGS *gs;
3143542a6bcSPeter Brune 
3153542a6bcSPeter Brune   PetscFunctionBegin;
316be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
317be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
318be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
319be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
320be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
321be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3223542a6bcSPeter Brune 
3233542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
324efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
3253542a6bcSPeter Brune 
3264fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
3274fc747eaSLawrence Mitchell 
32888976e71SPeter Brune   if (!snes->tolerancesset) {
3290e444f03SPeter Brune     snes->max_its   = 10000;
3300e444f03SPeter Brune     snes->max_funcs = 10000;
33188976e71SPeter Brune   }
3320e444f03SPeter Brune 
3334dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&gs));
3346cd56b8bSPeter Brune 
3356cd56b8bSPeter Brune   gs->sweeps  = 1;
336b6266c6eSPeter Brune   gs->rtol    = 1e-5;
33778e24b28SBarry Smith   gs->abstol  = PETSC_MACHINE_EPSILON;
33878e24b28SBarry Smith   gs->stol    = 1000 * PETSC_MACHINE_EPSILON;
339b6266c6eSPeter Brune   gs->max_its = 50;
34078e24b28SBarry Smith   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
3416cd56b8bSPeter Brune 
3423542a6bcSPeter Brune   snes->data = (void *)gs;
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3443542a6bcSPeter Brune }
345