xref: /petsc/src/snes/impls/gs/snesgs.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h"  I*/
26cd56b8bSPeter Brune 
3b6266c6eSPeter Brune /*@
4be95d8f1SBarry Smith    SNESNGSSetTolerances - Sets various parameters used in convergence tests.
5b6266c6eSPeter Brune 
6b6266c6eSPeter Brune    Logically Collective on SNES
7b6266c6eSPeter Brune 
8b6266c6eSPeter Brune    Input Parameters:
9b6266c6eSPeter Brune +  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 
23db781477SPatrick Sanan .seealso: `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 /*@
51be95d8f1SBarry Smith    SNESNGSGetTolerances - Gets various parameters used in convergence tests.
52b6266c6eSPeter Brune 
53b6266c6eSPeter Brune    Not Collective
54b6266c6eSPeter Brune 
55b6266c6eSPeter Brune    Input Parameters:
56b6266c6eSPeter Brune +  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 
63b6266c6eSPeter Brune    Notes:
640298fd71SBarry Smith    The user can specify NULL for any parameter that is not needed.
65b6266c6eSPeter Brune 
66b6266c6eSPeter Brune    Level: intermediate
67b6266c6eSPeter Brune 
68db781477SPatrick Sanan .seealso: `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 /*@
83be95d8f1SBarry Smith    SNESNGSSetSweeps - Sets the number of sweeps of GS to use.
846cd56b8bSPeter Brune 
856cd56b8bSPeter Brune    Input Parameters:
866cd56b8bSPeter Brune +  snes   - the SNES context
876cd56b8bSPeter Brune -  sweeps  - the number of sweeps of GS to perform.
886cd56b8bSPeter Brune 
896cd56b8bSPeter Brune    Level: intermediate
906cd56b8bSPeter Brune 
91db781477SPatrick Sanan .seealso: `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()`
926cd56b8bSPeter Brune @*/
936cd56b8bSPeter Brune 
949371c9d4SSatish Balay PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps) {
95be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
966cd56b8bSPeter Brune 
976cd56b8bSPeter Brune   PetscFunctionBegin;
986cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
996cd56b8bSPeter Brune   gs->sweeps = sweeps;
1006cd56b8bSPeter Brune   PetscFunctionReturn(0);
1016cd56b8bSPeter Brune }
1026cd56b8bSPeter Brune 
1036cd56b8bSPeter Brune /*@
104be95d8f1SBarry Smith    SNESNGSGetSweeps - Gets the number of sweeps GS will use.
1056cd56b8bSPeter Brune 
1066cd56b8bSPeter Brune    Input Parameters:
1076cd56b8bSPeter Brune .  snes   - the SNES context
1086cd56b8bSPeter Brune 
1096cd56b8bSPeter Brune    Output Parameters:
1106cd56b8bSPeter Brune .  sweeps  - the number of sweeps of GS to perform.
1116cd56b8bSPeter Brune 
1126cd56b8bSPeter Brune    Level: intermediate
1136cd56b8bSPeter Brune 
114db781477SPatrick Sanan .seealso: `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()`
1156cd56b8bSPeter Brune @*/
1169371c9d4SSatish Balay PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps) {
117be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1186cd56b8bSPeter Brune 
1196cd56b8bSPeter Brune   PetscFunctionBegin;
1206cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1216cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1226cd56b8bSPeter Brune   PetscFunctionReturn(0);
1236cd56b8bSPeter Brune }
1246cd56b8bSPeter Brune 
1259371c9d4SSatish Balay PetscErrorCode SNESReset_NGS(SNES snes) {
1268a86d3c5SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
1278a86d3c5SBarry Smith 
1283542a6bcSPeter Brune   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(ISColoringDestroy(&gs->coloring));
1303542a6bcSPeter Brune   PetscFunctionReturn(0);
1313542a6bcSPeter Brune }
1323542a6bcSPeter Brune 
1339371c9d4SSatish Balay PetscErrorCode SNESDestroy_NGS(SNES snes) {
1343542a6bcSPeter Brune   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(SNESReset_NGS(snes));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
1373542a6bcSPeter Brune   PetscFunctionReturn(0);
1383542a6bcSPeter Brune }
1393542a6bcSPeter Brune 
1409371c9d4SSatish Balay PetscErrorCode SNESSetUp_NGS(SNES snes) {
141ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
142ef107cf5SPeter Brune 
1433542a6bcSPeter Brune   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(SNESGetNGS(snes, &f, NULL));
145*48a46eb9SPierre Jolivet   if (!f) PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1463542a6bcSPeter Brune   PetscFunctionReturn(0);
1473542a6bcSPeter Brune }
1483542a6bcSPeter Brune 
1499371c9d4SSatish Balay PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems *PetscOptionsObject) {
150be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
151b6266c6eSPeter Brune   PetscInt  sweeps, max_its = PETSC_DEFAULT;
152b6266c6eSPeter Brune   PetscReal rtol = PETSC_DEFAULT, atol = PETSC_DEFAULT, stol = PETSC_DEFAULT;
153b6266c6eSPeter Brune   PetscBool flg, flg1, flg2, flg3;
1543542a6bcSPeter Brune 
1553542a6bcSPeter Brune   PetscFunctionBegin;
156d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options");
1576cd56b8bSPeter Brune   /* GS Options */
1589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg));
1591baa6e33SBarry Smith   if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps));
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg));
1619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1));
1629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2));
1639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3));
164*48a46eb9SPierre Jolivet   if (flg || flg1 || flg2 || flg3) PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its));
165ef107cf5SPeter Brune   flg = PETSC_FALSE;
1669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL));
167ef107cf5SPeter Brune   if (flg) {
1689566063dSJacob Faibussowitsch     PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
1699566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n"));
170ef107cf5SPeter Brune   }
1719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL));
1729566063dSJacob 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));
173ef107cf5SPeter Brune 
174d0609cedSBarry Smith   PetscOptionsHeadEnd();
1753542a6bcSPeter Brune   PetscFunctionReturn(0);
1763542a6bcSPeter Brune }
1773542a6bcSPeter Brune 
1789371c9d4SSatish Balay PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer) {
17998272f38SBarry Smith   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
18098272f38SBarry Smith   SNES_NGS *gs = (SNES_NGS *)snes->data;
18198272f38SBarry Smith   PetscBool iascii;
18298272f38SBarry Smith 
1833542a6bcSPeter Brune   PetscFunctionBegin;
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
18598272f38SBarry Smith   if (iascii) {
1869566063dSJacob Faibussowitsch     PetscCall(DMSNESGetNGS(snes->dm, &f, NULL));
187*48a46eb9SPierre Jolivet     if (f == SNESComputeNGSDefaultSecant) PetscCall(PetscViewerASCIIPrintf(viewer, "  Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h));
18898272f38SBarry Smith   }
1893542a6bcSPeter Brune   PetscFunctionReturn(0);
1903542a6bcSPeter Brune }
1913542a6bcSPeter Brune 
1929371c9d4SSatish Balay PetscErrorCode SNESSolve_NGS(SNES snes) {
1933542a6bcSPeter Brune   Vec              F;
1943542a6bcSPeter Brune   Vec              X;
1953542a6bcSPeter Brune   Vec              B;
1963542a6bcSPeter Brune   PetscInt         i;
19706e07b1aSPeter Brune   PetscReal        fnorm;
198365a6726SPeter Brune   SNESNormSchedule normschedule;
1993542a6bcSPeter Brune 
2003542a6bcSPeter Brune   PetscFunctionBegin;
201c579b300SPatrick Farrell 
2020b121fc5SBarry 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);
203c579b300SPatrick Farrell 
2049566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
2053542a6bcSPeter Brune   X = snes->vec_sol;
2063542a6bcSPeter Brune   F = snes->vec_func;
2073542a6bcSPeter Brune   B = snes->vec_rhs;
2083542a6bcSPeter Brune 
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
210534ebe21SPeter Brune   snes->iter = 0;
211534ebe21SPeter Brune   snes->norm = 0.;
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
213526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
214534ebe21SPeter Brune 
2159566063dSJacob Faibussowitsch   PetscCall(SNESGetNormSchedule(snes, &normschedule));
216365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2173542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
218e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2199566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2201aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
221e4ed7901SPeter Brune 
2229566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
223422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2249566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
225e4ed7901SPeter Brune     snes->iter = 0;
2263542a6bcSPeter Brune     snes->norm = fnorm;
2279566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2289566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
2299566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, 0, snes->norm));
2303542a6bcSPeter Brune 
2313542a6bcSPeter Brune     /* test convergence */
232dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
2333542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
234534ebe21SPeter Brune   } else {
2359566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2369566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
237534ebe21SPeter Brune   }
238534ebe21SPeter Brune 
2393542a6bcSPeter Brune   /* Call general purpose update function */
240dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
241534ebe21SPeter Brune 
2423542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
2439566063dSJacob Faibussowitsch     PetscCall(SNESComputeNGS(snes, B, X));
2444b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
245365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2469566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2479566063dSJacob Faibussowitsch       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
248422a814eSBarry Smith       SNESCheckFunctionNorm(snes, fnorm);
2493542a6bcSPeter Brune       /* Monitor convergence */
2509566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
2513542a6bcSPeter Brune       snes->iter = i + 1;
2523542a6bcSPeter Brune       snes->norm = fnorm;
2539566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2549566063dSJacob Faibussowitsch       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
2559566063dSJacob Faibussowitsch       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
256ac4a6415SMatthew G. Knepley     }
2573542a6bcSPeter Brune     /* Test for convergence */
258dbbe0bcdSBarry Smith     if (normschedule == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
2593542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2603542a6bcSPeter Brune     /* Call general purpose update function */
261dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
2623542a6bcSPeter Brune   }
263365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
264534ebe21SPeter Brune     if (i == snes->max_its) {
26563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its));
266534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
267534ebe21SPeter Brune     }
2681aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
2693542a6bcSPeter Brune   PetscFunctionReturn(0);
2703542a6bcSPeter Brune }
2713542a6bcSPeter Brune 
2723542a6bcSPeter Brune /*MC
27398272f38SBarry Smith   SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation
27498272f38SBarry Smith             using coloring.
2753542a6bcSPeter Brune 
2763542a6bcSPeter Brune    Level: advanced
2773542a6bcSPeter Brune 
27898272f38SBarry Smith   Options Database:
27978e24b28SBarry Smith +   -snes_ngs_sweeps <n> - Number of sweeps of GS to apply
28078e24b28SBarry Smith .    -snes_ngs_atol <atol> - Absolute residual tolerance for GS iteration
28178e24b28SBarry Smith .    -snes_ngs_rtol <rtol> - Relative residual tolerance for GS iteration
28278e24b28SBarry Smith .    -snes_ngs_stol <stol> - Absolute update tolerance for GS iteration
28378e24b28SBarry Smith .    -snes_ngs_max_it <maxit> - Maximum number of sweeps of GS to apply
28498272f38SBarry Smith .    -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is
28598272f38SBarry Smith                         used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring
28698272f38SBarry Smith                         is available or a Jacobian sparse matrix is provided (from which to get the coloring).
28778e24b28SBarry Smith .    -snes_ngs_secant_h <h> - Differencing parameter for secant approximation
28878e24b28SBarry Smith .    -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available.
289a5b23f4aSJose E. Roman -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed
29098272f38SBarry Smith 
2913542a6bcSPeter Brune   Notes:
292d814e60bSBarry Smith   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetNPC(), it will have
2933542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
2943542a6bcSPeter Brune 
29578e24b28SBarry Smith   By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with SNESSetNormSchedule() or -snes_norm_schedule
2964f02bc6aSBarry Smith    References:
297606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
2984f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
2994f02bc6aSBarry Smith 
300db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`,
301db781477SPatrick Sanan           `SNESSetNormSchedule()`
3023542a6bcSPeter Brune M*/
3033542a6bcSPeter Brune 
3049371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes) {
305be95d8f1SBarry Smith   SNES_NGS *gs;
3063542a6bcSPeter Brune 
3073542a6bcSPeter Brune   PetscFunctionBegin;
308be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
309be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
310be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
311be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
312be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
313be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3143542a6bcSPeter Brune 
3153542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
316efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
3173542a6bcSPeter Brune 
3184fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
3194fc747eaSLawrence Mitchell 
32088976e71SPeter Brune   if (!snes->tolerancesset) {
3210e444f03SPeter Brune     snes->max_its   = 10000;
3220e444f03SPeter Brune     snes->max_funcs = 10000;
32388976e71SPeter Brune   }
3240e444f03SPeter Brune 
3259566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes, &gs));
3266cd56b8bSPeter Brune 
3276cd56b8bSPeter Brune   gs->sweeps  = 1;
328b6266c6eSPeter Brune   gs->rtol    = 1e-5;
32978e24b28SBarry Smith   gs->abstol  = PETSC_MACHINE_EPSILON;
33078e24b28SBarry Smith   gs->stol    = 1000 * PETSC_MACHINE_EPSILON;
331b6266c6eSPeter Brune   gs->max_its = 50;
33278e24b28SBarry Smith   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
3336cd56b8bSPeter Brune 
3343542a6bcSPeter Brune   snes->data = (void *)gs;
3353542a6bcSPeter Brune   PetscFunctionReturn(0);
3363542a6bcSPeter Brune }
337