xref: /petsc/src/snes/impls/gs/snesgs.c (revision f6dfbefd03961ab3be6f06be75c96cbf27a49667)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h"  I*/
26cd56b8bSPeter Brune 
3b6266c6eSPeter Brune /*@
4*f6dfbefdSBarry Smith    SNESNGSSetTolerances - Sets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG`
5b6266c6eSPeter Brune 
6*f6dfbefdSBarry Smith    Logically Collective on snes
7b6266c6eSPeter Brune 
8b6266c6eSPeter Brune    Input Parameters:
9*f6dfbefdSBarry 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 
23*f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetTrustRegionTolerance()`
24b6266c6eSPeter Brune @*/
259371c9d4SSatish Balay PetscErrorCode SNESNGSSetTolerances(SNES snes, PetscReal abstol, PetscReal rtol, PetscReal stol, PetscInt maxit) {
26be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
27b6266c6eSPeter Brune 
28b6266c6eSPeter Brune   PetscFunctionBegin;
29b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
30b6266c6eSPeter Brune 
31b6266c6eSPeter Brune   if (abstol != PETSC_DEFAULT) {
3208401ef6SPierre Jolivet     PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
33b6266c6eSPeter Brune     gs->abstol = abstol;
34b6266c6eSPeter Brune   }
35b6266c6eSPeter Brune   if (rtol != PETSC_DEFAULT) {
360b121fc5SBarry 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);
37b6266c6eSPeter Brune     gs->rtol = rtol;
38b6266c6eSPeter Brune   }
39b6266c6eSPeter Brune   if (stol != PETSC_DEFAULT) {
4008401ef6SPierre Jolivet     PetscCheck(stol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Step tolerance %g must be non-negative", (double)stol);
41b6266c6eSPeter Brune     gs->stol = stol;
42b6266c6eSPeter Brune   }
43b6266c6eSPeter Brune   if (maxit != PETSC_DEFAULT) {
4463a3b9bcSJacob Faibussowitsch     PetscCheck(maxit >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxit);
45b6266c6eSPeter Brune     gs->max_its = maxit;
46b6266c6eSPeter Brune   }
47b6266c6eSPeter Brune   PetscFunctionReturn(0);
48b6266c6eSPeter Brune }
49b6266c6eSPeter Brune 
50b6266c6eSPeter Brune /*@
51*f6dfbefdSBarry Smith    SNESNGSGetTolerances - Gets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG`
52b6266c6eSPeter Brune 
53b6266c6eSPeter Brune    Not Collective
54b6266c6eSPeter Brune 
55b6266c6eSPeter Brune    Input Parameters:
56*f6dfbefdSBarry Smith +  snes - the `SNES` context
57b6266c6eSPeter Brune .  atol - absolute convergence tolerance
58b6266c6eSPeter Brune .  rtol - relative convergence tolerance
59b6266c6eSPeter Brune .  stol -  convergence tolerance in terms of the norm
60b6266c6eSPeter Brune            of the change in the solution between steps
61b6266c6eSPeter Brune -  maxit - maximum number of iterations
62b6266c6eSPeter Brune 
63*f6dfbefdSBarry Smith    Note:
640298fd71SBarry Smith    The user can specify NULL for any parameter that is not needed.
65b6266c6eSPeter Brune 
66b6266c6eSPeter Brune    Level: intermediate
67b6266c6eSPeter Brune 
68*f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetTolerances()`
69b6266c6eSPeter Brune @*/
709371c9d4SSatish Balay PetscErrorCode SNESNGSGetTolerances(SNES snes, PetscReal *atol, PetscReal *rtol, PetscReal *stol, PetscInt *maxit) {
71be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
72b6266c6eSPeter Brune 
73b6266c6eSPeter Brune   PetscFunctionBegin;
74b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
75b6266c6eSPeter Brune   if (atol) *atol = gs->abstol;
76b6266c6eSPeter Brune   if (rtol) *rtol = gs->rtol;
77b6266c6eSPeter Brune   if (stol) *stol = gs->stol;
78b6266c6eSPeter Brune   if (maxit) *maxit = gs->max_its;
79b6266c6eSPeter Brune   PetscFunctionReturn(0);
80b6266c6eSPeter Brune }
816cd56b8bSPeter Brune 
826cd56b8bSPeter Brune /*@
83*f6dfbefdSBarry Smith    SNESNGSSetSweeps - Sets the number of sweeps of nonlinear GS to use in `SNESNCG`
846cd56b8bSPeter Brune 
856cd56b8bSPeter Brune    Input Parameters:
86*f6dfbefdSBarry Smith +  snes   - the `SNES` context
87*f6dfbefdSBarry Smith -  sweeps  - the number of sweeps of nonlinear GS to perform.
88*f6dfbefdSBarry Smith 
89*f6dfbefdSBarry Smith   Options Database Key:
90*f6dfbefdSBarry Smith .   -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply
916cd56b8bSPeter Brune 
926cd56b8bSPeter Brune    Level: intermediate
936cd56b8bSPeter Brune 
94*f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()`
956cd56b8bSPeter Brune @*/
966cd56b8bSPeter Brune 
979371c9d4SSatish Balay PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps) {
98be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
996cd56b8bSPeter Brune 
1006cd56b8bSPeter Brune   PetscFunctionBegin;
1016cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1026cd56b8bSPeter Brune   gs->sweeps = sweeps;
1036cd56b8bSPeter Brune   PetscFunctionReturn(0);
1046cd56b8bSPeter Brune }
1056cd56b8bSPeter Brune 
1066cd56b8bSPeter Brune /*@
107*f6dfbefdSBarry Smith    SNESNGSGetSweeps - Gets the number of sweeps nonlinear GS will use in `SNESNCG`
1086cd56b8bSPeter Brune 
1096cd56b8bSPeter Brune    Input Parameters:
110*f6dfbefdSBarry Smith .  snes   - the `SNES` context
1116cd56b8bSPeter Brune 
1126cd56b8bSPeter Brune    Output Parameters:
113*f6dfbefdSBarry Smith .  sweeps  - the number of sweeps of nonlinear GS to perform.
1146cd56b8bSPeter Brune 
1156cd56b8bSPeter Brune    Level: intermediate
1166cd56b8bSPeter Brune 
117*f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()`
1186cd56b8bSPeter Brune @*/
1199371c9d4SSatish Balay PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps) {
120be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1216cd56b8bSPeter Brune 
1226cd56b8bSPeter Brune   PetscFunctionBegin;
1236cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1246cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1256cd56b8bSPeter Brune   PetscFunctionReturn(0);
1266cd56b8bSPeter Brune }
1276cd56b8bSPeter Brune 
1289371c9d4SSatish Balay PetscErrorCode SNESReset_NGS(SNES snes) {
1298a86d3c5SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1308a86d3c5SBarry Smith 
1313542a6bcSPeter Brune   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCall(ISColoringDestroy(&gs->coloring));
1333542a6bcSPeter Brune   PetscFunctionReturn(0);
1343542a6bcSPeter Brune }
1353542a6bcSPeter Brune 
1369371c9d4SSatish Balay PetscErrorCode SNESDestroy_NGS(SNES snes) {
1373542a6bcSPeter Brune   PetscFunctionBegin;
1389566063dSJacob Faibussowitsch   PetscCall(SNESReset_NGS(snes));
1399566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
1403542a6bcSPeter Brune   PetscFunctionReturn(0);
1413542a6bcSPeter Brune }
1423542a6bcSPeter Brune 
1439371c9d4SSatish Balay PetscErrorCode SNESSetUp_NGS(SNES snes) {
144ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
145ef107cf5SPeter Brune 
1463542a6bcSPeter Brune   PetscFunctionBegin;
1479566063dSJacob Faibussowitsch   PetscCall(SNESGetNGS(snes, &f, NULL));
14848a46eb9SPierre Jolivet   if (!f) PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1493542a6bcSPeter Brune   PetscFunctionReturn(0);
1503542a6bcSPeter Brune }
1513542a6bcSPeter Brune 
1529371c9d4SSatish Balay PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems *PetscOptionsObject) {
153be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
154b6266c6eSPeter Brune   PetscInt  sweeps, max_its = PETSC_DEFAULT;
155b6266c6eSPeter Brune   PetscReal rtol = PETSC_DEFAULT, atol = PETSC_DEFAULT, stol = PETSC_DEFAULT;
156b6266c6eSPeter Brune   PetscBool flg, flg1, flg2, flg3;
1573542a6bcSPeter Brune 
1583542a6bcSPeter Brune   PetscFunctionBegin;
159d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options");
1606cd56b8bSPeter Brune   /* GS Options */
1619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg));
1621baa6e33SBarry Smith   if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps));
1639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg));
1649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1));
1659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2));
1669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3));
16748a46eb9SPierre Jolivet   if (flg || flg1 || flg2 || flg3) PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its));
168ef107cf5SPeter Brune   flg = PETSC_FALSE;
1699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL));
170ef107cf5SPeter Brune   if (flg) {
1719566063dSJacob Faibussowitsch     PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1729566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n"));
173ef107cf5SPeter Brune   }
1749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL));
1759566063dSJacob 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));
176ef107cf5SPeter Brune 
177d0609cedSBarry Smith   PetscOptionsHeadEnd();
1783542a6bcSPeter Brune   PetscFunctionReturn(0);
1793542a6bcSPeter Brune }
1803542a6bcSPeter Brune 
1819371c9d4SSatish Balay PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer) {
18298272f38SBarry Smith   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
18398272f38SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
18498272f38SBarry Smith   PetscBool iascii;
18598272f38SBarry Smith 
1863542a6bcSPeter Brune   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
18898272f38SBarry Smith   if (iascii) {
1899566063dSJacob Faibussowitsch     PetscCall(DMSNESGetNGS(snes->dm, &f, NULL));
19048a46eb9SPierre Jolivet     if (f == SNESComputeNGSDefaultSecant) PetscCall(PetscViewerASCIIPrintf(viewer, "  Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h));
19198272f38SBarry Smith   }
1923542a6bcSPeter Brune   PetscFunctionReturn(0);
1933542a6bcSPeter Brune }
1943542a6bcSPeter Brune 
1959371c9d4SSatish Balay PetscErrorCode SNESSolve_NGS(SNES snes) {
1963542a6bcSPeter Brune   Vec              F;
1973542a6bcSPeter Brune   Vec              X;
1983542a6bcSPeter Brune   Vec              B;
1993542a6bcSPeter Brune   PetscInt         i;
20006e07b1aSPeter Brune   PetscReal        fnorm;
201365a6726SPeter Brune   SNESNormSchedule normschedule;
2023542a6bcSPeter Brune 
2033542a6bcSPeter Brune   PetscFunctionBegin;
204c579b300SPatrick Farrell 
2050b121fc5SBarry 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);
206c579b300SPatrick Farrell 
2079566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
2083542a6bcSPeter Brune   X = snes->vec_sol;
2093542a6bcSPeter Brune   F = snes->vec_func;
2103542a6bcSPeter Brune   B = snes->vec_rhs;
2113542a6bcSPeter Brune 
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
213534ebe21SPeter Brune   snes->iter = 0;
214534ebe21SPeter Brune   snes->norm = 0.;
2159566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
216526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
217534ebe21SPeter Brune 
2189566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
219365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2203542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
221e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2229566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2231aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
224e4ed7901SPeter Brune 
2259566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
226422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2279566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
228e4ed7901SPeter Brune     snes->iter = 0;
2293542a6bcSPeter Brune     snes->norm = fnorm;
2309566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2319566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
2329566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, 0, snes->norm));
2333542a6bcSPeter Brune 
2343542a6bcSPeter Brune     /* test convergence */
235dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
2363542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
237534ebe21SPeter Brune   } else {
2389566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2399566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
240534ebe21SPeter Brune   }
241534ebe21SPeter Brune 
2423542a6bcSPeter Brune   /* Call general purpose update function */
243dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
244534ebe21SPeter Brune 
2453542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
2469566063dSJacob Faibussowitsch     PetscCall(SNESComputeNGS(snes, B, X));
2474b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
248365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2499566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2509566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
251422a814eSBarry Smith       SNESCheckFunctionNorm(snes, fnorm);
2523542a6bcSPeter Brune       /* Monitor convergence */
2539566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
2543542a6bcSPeter Brune       snes->iter = i + 1;
2553542a6bcSPeter Brune       snes->norm = fnorm;
2569566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2579566063dSJacob Faibussowitsch       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
2589566063dSJacob Faibussowitsch       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
259ac4a6415SMatthew G. Knepley     }
2603542a6bcSPeter Brune     /* Test for convergence */
261dbbe0bcdSBarry Smith     if (normschedule == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
2623542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2633542a6bcSPeter Brune     /* Call general purpose update function */
264dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
2653542a6bcSPeter Brune   }
266365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
267534ebe21SPeter Brune     if (i == snes->max_its) {
26863a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its));
269534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
270534ebe21SPeter Brune     }
2711aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
2723542a6bcSPeter Brune   PetscFunctionReturn(0);
2733542a6bcSPeter Brune }
2743542a6bcSPeter Brune 
2753542a6bcSPeter Brune /*MC
276*f6dfbefdSBarry Smith   SNESNGS - Either calls the user-provided solution routine provided with `SNESSetNGS()` or does a finite difference secant approximation
27798272f38SBarry Smith             using coloring.
2783542a6bcSPeter Brune 
2793542a6bcSPeter Brune    Level: advanced
2803542a6bcSPeter Brune 
281*f6dfbefdSBarry Smith   Options Database Keys:
282*f6dfbefdSBarry Smith +   -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply
283*f6dfbefdSBarry Smith .    -snes_ngs_atol <atol> - Absolute residual tolerance for nonlinear GS iteration
284*f6dfbefdSBarry Smith .    -snes_ngs_rtol <rtol> - Relative residual tolerance for nonlinear GS iteration
285*f6dfbefdSBarry Smith .    -snes_ngs_stol <stol> - Absolute update tolerance for nonlinear GS iteration
286*f6dfbefdSBarry Smith .    -snes_ngs_max_it <maxit> - Maximum number of sweeps of nonlinea GS to apply
28798272f38SBarry Smith .    -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is
288*f6dfbefdSBarry Smith                         used by default if no user provided Gauss-Seidel routine is available. Requires either that a `DM` that can compute a coloring
28998272f38SBarry Smith                         is available or a Jacobian sparse matrix is provided (from which to get the coloring).
29078e24b28SBarry Smith .    -snes_ngs_secant_h <h> - Differencing parameter for secant approximation
29178e24b28SBarry Smith .    -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available.
292a5b23f4aSJose E. Roman -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed
29398272f38SBarry Smith 
2943542a6bcSPeter Brune   Notes:
295*f6dfbefdSBarry Smith   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with `SNESGetNPC()`, it will have
2963542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
2973542a6bcSPeter Brune 
298*f6dfbefdSBarry Smith   By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with `SNESSetNormSchedule()`
299*f6dfbefdSBarry Smith   or -snes_norm_schedule none
300*f6dfbefdSBarry Smith 
3014f02bc6aSBarry Smith    References:
302606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
3034f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
3044f02bc6aSBarry Smith 
305*f6dfbefdSBarry Smith .seealso: `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`,
306db781477SPatrick Sanan           `SNESSetNormSchedule()`
3073542a6bcSPeter Brune M*/
3083542a6bcSPeter Brune 
3099371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes) {
310be95d8f1SBarry Smith   SNES_NGS *gs;
3113542a6bcSPeter Brune 
3123542a6bcSPeter Brune   PetscFunctionBegin;
313be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
314be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
315be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
316be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
317be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
318be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3193542a6bcSPeter Brune 
3203542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
321efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
3223542a6bcSPeter Brune 
3234fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
3244fc747eaSLawrence Mitchell 
32588976e71SPeter Brune   if (!snes->tolerancesset) {
3260e444f03SPeter Brune     snes->max_its   = 10000;
3270e444f03SPeter Brune     snes->max_funcs = 10000;
32888976e71SPeter Brune   }
3290e444f03SPeter Brune 
3309566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &gs));
3316cd56b8bSPeter Brune 
3326cd56b8bSPeter Brune   gs->sweeps  = 1;
333b6266c6eSPeter Brune   gs->rtol    = 1e-5;
33478e24b28SBarry Smith   gs->abstol  = PETSC_MACHINE_EPSILON;
33578e24b28SBarry Smith   gs->stol    = 1000 * PETSC_MACHINE_EPSILON;
336b6266c6eSPeter Brune   gs->max_its = 50;
33778e24b28SBarry Smith   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
3386cd56b8bSPeter Brune 
3393542a6bcSPeter Brune   snes->data = (void *)gs;
3403542a6bcSPeter Brune   PetscFunctionReturn(0);
3413542a6bcSPeter Brune }
342