xref: /petsc/src/snes/impls/gs/snesgs.c (revision 78e24b289582e8ff394a680fe2cef1e75b7ea512)
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 
16b6266c6eSPeter Brune    Options Database Keys:
17be95d8f1SBarry Smith +    -snes_ngs_atol <abstol> - Sets abstol
18be95d8f1SBarry Smith .    -snes_ngs_rtol <rtol> - Sets rtol
19be95d8f1SBarry Smith .    -snes_ngs_stol <stol> - Sets stol
20b6266c6eSPeter Brune -    -snes_max_it <maxit> - Sets maxit
21b6266c6eSPeter Brune 
22b6266c6eSPeter Brune    Level: intermediate
23b6266c6eSPeter Brune 
24b6266c6eSPeter Brune .seealso: SNESSetTrustRegionTolerance()
25b6266c6eSPeter Brune @*/
26be95d8f1SBarry Smith PetscErrorCode  SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
27b6266c6eSPeter Brune {
28be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
29b6266c6eSPeter Brune 
30b6266c6eSPeter Brune   PetscFunctionBegin;
31b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
32b6266c6eSPeter Brune 
33b6266c6eSPeter Brune   if (abstol != PETSC_DEFAULT) {
3457622a8eSBarry Smith     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
35b6266c6eSPeter Brune     gs->abstol = abstol;
36b6266c6eSPeter Brune   }
37b6266c6eSPeter Brune   if (rtol != PETSC_DEFAULT) {
3857622a8eSBarry Smith     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol);
39b6266c6eSPeter Brune     gs->rtol = rtol;
40b6266c6eSPeter Brune   }
41b6266c6eSPeter Brune   if (stol != PETSC_DEFAULT) {
4257622a8eSBarry Smith     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
43b6266c6eSPeter Brune     gs->stol = stol;
44b6266c6eSPeter Brune   }
45b6266c6eSPeter Brune   if (maxit != PETSC_DEFAULT) {
46ce94432eSBarry Smith     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
47b6266c6eSPeter Brune     gs->max_its = maxit;
48b6266c6eSPeter Brune   }
49b6266c6eSPeter Brune   PetscFunctionReturn(0);
50b6266c6eSPeter Brune }
51b6266c6eSPeter Brune 
52b6266c6eSPeter Brune /*@
53be95d8f1SBarry Smith    SNESNGSGetTolerances - Gets various parameters used in convergence tests.
54b6266c6eSPeter Brune 
55b6266c6eSPeter Brune    Not Collective
56b6266c6eSPeter Brune 
57b6266c6eSPeter Brune    Input Parameters:
58b6266c6eSPeter Brune +  snes - the SNES context
59b6266c6eSPeter Brune .  atol - absolute convergence tolerance
60b6266c6eSPeter Brune .  rtol - relative convergence tolerance
61b6266c6eSPeter Brune .  stol -  convergence tolerance in terms of the norm
62b6266c6eSPeter Brune            of the change in the solution between steps
63b6266c6eSPeter Brune -  maxit - maximum number of iterations
64b6266c6eSPeter Brune 
65b6266c6eSPeter Brune    Notes:
660298fd71SBarry Smith    The user can specify NULL for any parameter that is not needed.
67b6266c6eSPeter Brune 
68b6266c6eSPeter Brune    Level: intermediate
69b6266c6eSPeter Brune 
70b6266c6eSPeter Brune .seealso: SNESSetTolerances()
71b6266c6eSPeter Brune @*/
72be95d8f1SBarry Smith PetscErrorCode  SNESNGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
73b6266c6eSPeter Brune {
74be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
75b6266c6eSPeter Brune 
76b6266c6eSPeter Brune   PetscFunctionBegin;
77b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
78b6266c6eSPeter Brune   if (atol) *atol = gs->abstol;
79b6266c6eSPeter Brune   if (rtol) *rtol = gs->rtol;
80b6266c6eSPeter Brune   if (stol) *stol = gs->stol;
81b6266c6eSPeter Brune   if (maxit) *maxit = gs->max_its;
82b6266c6eSPeter Brune   PetscFunctionReturn(0);
83b6266c6eSPeter Brune }
846cd56b8bSPeter Brune 
856cd56b8bSPeter Brune /*@
86be95d8f1SBarry Smith    SNESNGSSetSweeps - Sets the number of sweeps of GS to use.
876cd56b8bSPeter Brune 
886cd56b8bSPeter Brune    Input Parameters:
896cd56b8bSPeter Brune +  snes   - the SNES context
906cd56b8bSPeter Brune -  sweeps  - the number of sweeps of GS to perform.
916cd56b8bSPeter Brune 
926cd56b8bSPeter Brune    Level: intermediate
936cd56b8bSPeter Brune 
94d814e60bSBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetNPC(), SNESNGSGetSweeps()
956cd56b8bSPeter Brune @*/
966cd56b8bSPeter Brune 
97be95d8f1SBarry Smith PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
986cd56b8bSPeter Brune {
99be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
1006cd56b8bSPeter Brune 
1016cd56b8bSPeter Brune   PetscFunctionBegin;
1026cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1036cd56b8bSPeter Brune   gs->sweeps = sweeps;
1046cd56b8bSPeter Brune   PetscFunctionReturn(0);
1056cd56b8bSPeter Brune }
1066cd56b8bSPeter Brune 
1076cd56b8bSPeter Brune /*@
108be95d8f1SBarry Smith    SNESNGSGetSweeps - Gets the number of sweeps GS will use.
1096cd56b8bSPeter Brune 
1106cd56b8bSPeter Brune    Input Parameters:
1116cd56b8bSPeter Brune .  snes   - the SNES context
1126cd56b8bSPeter Brune 
1136cd56b8bSPeter Brune    Output Parameters:
1146cd56b8bSPeter Brune .  sweeps  - the number of sweeps of GS to perform.
1156cd56b8bSPeter Brune 
1166cd56b8bSPeter Brune    Level: intermediate
1176cd56b8bSPeter Brune 
118d814e60bSBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetNPC(), SNESNGSSetSweeps()
1196cd56b8bSPeter Brune @*/
120be95d8f1SBarry Smith PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt * sweeps)
1216cd56b8bSPeter Brune {
122be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
1236cd56b8bSPeter Brune 
1246cd56b8bSPeter Brune   PetscFunctionBegin;
1256cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1266cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1276cd56b8bSPeter Brune   PetscFunctionReturn(0);
1286cd56b8bSPeter Brune }
1296cd56b8bSPeter Brune 
130be95d8f1SBarry Smith PetscErrorCode SNESReset_NGS(SNES snes)
1313542a6bcSPeter Brune {
1328a86d3c5SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
1338a86d3c5SBarry Smith   PetscErrorCode ierr;
1348a86d3c5SBarry Smith 
1353542a6bcSPeter Brune   PetscFunctionBegin;
1368a86d3c5SBarry Smith   ierr = ISColoringDestroy(&gs->coloring);CHKERRQ(ierr);
1373542a6bcSPeter Brune   PetscFunctionReturn(0);
1383542a6bcSPeter Brune }
1393542a6bcSPeter Brune 
140be95d8f1SBarry Smith PetscErrorCode SNESDestroy_NGS(SNES snes)
1413542a6bcSPeter Brune {
1423542a6bcSPeter Brune   PetscErrorCode ierr;
1433542a6bcSPeter Brune 
1443542a6bcSPeter Brune   PetscFunctionBegin;
145be95d8f1SBarry Smith   ierr = SNESReset_NGS(snes);CHKERRQ(ierr);
14622d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1473542a6bcSPeter Brune   PetscFunctionReturn(0);
1483542a6bcSPeter Brune }
1493542a6bcSPeter Brune 
150be95d8f1SBarry Smith PetscErrorCode SNESSetUp_NGS(SNES snes)
1513542a6bcSPeter Brune {
152ef107cf5SPeter Brune   PetscErrorCode ierr;
153ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
154ef107cf5SPeter Brune 
1553542a6bcSPeter Brune   PetscFunctionBegin;
156be95d8f1SBarry Smith   ierr = SNESGetNGS(snes,&f,NULL);CHKERRQ(ierr);
157ef107cf5SPeter Brune   if (!f) {
158be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
159ef107cf5SPeter Brune   }
1603542a6bcSPeter Brune   PetscFunctionReturn(0);
1613542a6bcSPeter Brune }
1623542a6bcSPeter Brune 
1634416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGS(PetscOptionItems *PetscOptionsObject,SNES snes)
1643542a6bcSPeter Brune {
165be95d8f1SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
1663542a6bcSPeter Brune   PetscErrorCode ierr;
167b6266c6eSPeter Brune   PetscInt       sweeps,max_its=PETSC_DEFAULT;
168b6266c6eSPeter Brune   PetscReal      rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT;
169b6266c6eSPeter Brune   PetscBool      flg,flg1,flg2,flg3;
1703542a6bcSPeter Brune 
1713542a6bcSPeter Brune   PetscFunctionBegin;
172e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES GS options");CHKERRQ(ierr);
1736cd56b8bSPeter Brune   /* GS Options */
174be95d8f1SBarry Smith   ierr = PetscOptionsInt("-snes_ngs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr);
175b6266c6eSPeter Brune   if (flg) {
176be95d8f1SBarry Smith     ierr = SNESNGSSetSweeps(snes,sweeps);CHKERRQ(ierr);
177b6266c6eSPeter Brune   }
17851b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr);
17951b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr);
18051b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr);
18151b6f3e6SMatthew G. Knepley   ierr = PetscOptionsInt("-snes_ngs_max_it","Maximum number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);CHKERRQ(ierr);
182b6266c6eSPeter Brune   if (flg || flg1 || flg2 || flg3) {
183be95d8f1SBarry Smith     ierr = SNESNGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr);
184b6266c6eSPeter Brune   }
185ef107cf5SPeter Brune   flg  = PETSC_FALSE;
18698272f38SBarry Smith   ierr = PetscOptionsBool("-snes_ngs_secant","Use finite difference secant approximation with coloring","",flg,&flg,NULL);CHKERRQ(ierr);
187ef107cf5SPeter Brune   if (flg) {
188be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
18998272f38SBarry Smith     ierr = PetscInfo(snes,"Setting default finite difference secant approximation with coloring\n");CHKERRQ(ierr);
190ef107cf5SPeter Brune   }
191be95d8f1SBarry Smith   ierr = PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);CHKERRQ(ierr);
19298272f38SBarry Smith   ierr = PetscOptionsBool("-snes_ngs_secant_mat_coloring","Use the graph coloring of the Jacobian for the secant GS","",gs->secant_mat,&gs->secant_mat,&flg);CHKERRQ(ierr);
193ef107cf5SPeter Brune 
1944b32a720SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
1953542a6bcSPeter Brune   PetscFunctionReturn(0);
1963542a6bcSPeter Brune }
1973542a6bcSPeter Brune 
198be95d8f1SBarry Smith PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
1993542a6bcSPeter Brune {
20098272f38SBarry Smith   PetscErrorCode ierr;
20198272f38SBarry Smith   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
20298272f38SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
20398272f38SBarry Smith   PetscBool      iascii;
20498272f38SBarry Smith 
2053542a6bcSPeter Brune   PetscFunctionBegin;
20698272f38SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
20798272f38SBarry Smith   if (iascii) {
20898272f38SBarry Smith     ierr = DMSNESGetNGS(snes->dm,&f,NULL);CHKERRQ(ierr);
20998272f38SBarry Smith     if (f == SNESComputeNGSDefaultSecant) {
210efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Use finite difference secant approximation with coloring with h = %g \n",(double)gs->h);CHKERRQ(ierr);
21198272f38SBarry Smith     }
21298272f38SBarry Smith   }
2133542a6bcSPeter Brune   PetscFunctionReturn(0);
2143542a6bcSPeter Brune }
2153542a6bcSPeter Brune 
216be95d8f1SBarry Smith PetscErrorCode SNESSolve_NGS(SNES snes)
2173542a6bcSPeter Brune {
2183542a6bcSPeter Brune   Vec              F;
2193542a6bcSPeter Brune   Vec              X;
2203542a6bcSPeter Brune   Vec              B;
2213542a6bcSPeter Brune   PetscInt         i;
22206e07b1aSPeter Brune   PetscReal        fnorm;
2233542a6bcSPeter Brune   PetscErrorCode   ierr;
224365a6726SPeter Brune   SNESNormSchedule normschedule;
2253542a6bcSPeter Brune 
2263542a6bcSPeter Brune   PetscFunctionBegin;
227c579b300SPatrick Farrell 
2286c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
229c579b300SPatrick Farrell 
230fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
2313542a6bcSPeter Brune   X = snes->vec_sol;
2323542a6bcSPeter Brune   F = snes->vec_func;
2333542a6bcSPeter Brune   B = snes->vec_rhs;
2343542a6bcSPeter Brune 
235e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
236534ebe21SPeter Brune   snes->iter   = 0;
237534ebe21SPeter Brune   snes->norm   = 0.;
238e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
239526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
240534ebe21SPeter Brune 
241365a6726SPeter Brune   ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
242365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2433542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
244e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2453542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2461aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
247e4ed7901SPeter Brune 
2483542a6bcSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
249422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
250e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
251e4ed7901SPeter Brune     snes->iter = 0;
2523542a6bcSPeter Brune     snes->norm = fnorm;
253e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
254a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
255e4ed7901SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
2563542a6bcSPeter Brune 
2573542a6bcSPeter Brune     /* test convergence */
2583542a6bcSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2593542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
260534ebe21SPeter Brune   } else {
261e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
262a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
263534ebe21SPeter Brune   }
264534ebe21SPeter Brune 
2653542a6bcSPeter Brune   /* Call general purpose update function */
2663542a6bcSPeter Brune   if (snes->ops->update) {
2673542a6bcSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2683542a6bcSPeter Brune   }
269534ebe21SPeter Brune 
2703542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
271be95d8f1SBarry Smith     ierr = SNESComputeNGS(snes, B, X);CHKERRQ(ierr);
2724b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
273365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2743542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2753542a6bcSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
276422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
2773542a6bcSPeter Brune       /* Monitor convergence */
278e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
2793542a6bcSPeter Brune       snes->iter = i+1;
2803542a6bcSPeter Brune       snes->norm = fnorm;
281e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
282a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
2833542a6bcSPeter Brune       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
284ac4a6415SMatthew G. Knepley     }
2853542a6bcSPeter Brune     /* Test for convergence */
286365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2873542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2883542a6bcSPeter Brune     /* Call general purpose update function */
2893542a6bcSPeter Brune     if (snes->ops->update) {
2903542a6bcSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2913542a6bcSPeter Brune     }
2923542a6bcSPeter Brune   }
293365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
294534ebe21SPeter Brune     if (i == snes->max_its) {
295534ebe21SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
296534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
297534ebe21SPeter Brune     }
2981aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
2993542a6bcSPeter Brune   PetscFunctionReturn(0);
3003542a6bcSPeter Brune }
3013542a6bcSPeter Brune 
3023542a6bcSPeter Brune /*MC
30398272f38SBarry Smith   SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation
30498272f38SBarry Smith             using coloring.
3053542a6bcSPeter Brune 
3063542a6bcSPeter Brune    Level: advanced
3073542a6bcSPeter Brune 
30898272f38SBarry Smith   Options Database:
309*78e24b28SBarry Smith +   -snes_ngs_sweeps <n> - Number of sweeps of GS to apply
310*78e24b28SBarry Smith .    -snes_ngs_atol <atol> - Absolute residual tolerance for GS iteration
311*78e24b28SBarry Smith .    -snes_ngs_rtol <rtol> - Relative residual tolerance for GS iteration
312*78e24b28SBarry Smith .    -snes_ngs_stol <stol> - Absolute update tolerance for GS iteration
313*78e24b28SBarry Smith .    -snes_ngs_max_it <maxit> - Maximum number of sweeps of GS to apply
31498272f38SBarry Smith .    -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is
31598272f38SBarry Smith                         used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring
31698272f38SBarry Smith                         is available or a Jacobian sparse matrix is provided (from which to get the coloring).
317*78e24b28SBarry Smith .    -snes_ngs_secant_h <h> - Differencing parameter for secant approximation
318*78e24b28SBarry Smith .    -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available.
319*78e24b28SBarry Smith -    -snes_norm_schedule <none, always, initialonly, finalonly, initalfinalonly> - how often the residual norms are computed
32098272f38SBarry Smith 
3213542a6bcSPeter Brune   Notes:
322d814e60bSBarry Smith   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetNPC(), it will have
3233542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
3243542a6bcSPeter Brune 
325*78e24b28SBarry 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
3264f02bc6aSBarry Smith    References:
32796a0c994SBarry Smith .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
3284f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
3294f02bc6aSBarry Smith 
330*78e24b28SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetNGS(), SNESType (for list of available types), SNESNGSSetSweeps(), SNESNGSSetTolerances(),
331*78e24b28SBarry Smith           SNESSetNormSchedule()
3323542a6bcSPeter Brune M*/
3333542a6bcSPeter Brune 
334be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
3353542a6bcSPeter Brune {
336be95d8f1SBarry Smith   SNES_NGS        *gs;
3373542a6bcSPeter Brune   PetscErrorCode ierr;
3383542a6bcSPeter Brune 
3393542a6bcSPeter Brune   PetscFunctionBegin;
340be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
341be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
342be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
343be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
344be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
345be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3463542a6bcSPeter Brune 
3473542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
348efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
3493542a6bcSPeter Brune 
3504fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
3514fc747eaSLawrence Mitchell 
35288976e71SPeter Brune   if (!snes->tolerancesset) {
3530e444f03SPeter Brune     snes->max_its   = 10000;
3540e444f03SPeter Brune     snes->max_funcs = 10000;
35588976e71SPeter Brune   }
3560e444f03SPeter Brune 
357b00a9115SJed Brown   ierr = PetscNewLog(snes,&gs);CHKERRQ(ierr);
3586cd56b8bSPeter Brune 
3596cd56b8bSPeter Brune   gs->sweeps  = 1;
360b6266c6eSPeter Brune   gs->rtol    = 1e-5;
361*78e24b28SBarry Smith   gs->abstol  = PETSC_MACHINE_EPSILON;
362*78e24b28SBarry Smith   gs->stol    = 1000*PETSC_MACHINE_EPSILON;
363b6266c6eSPeter Brune   gs->max_its = 50;
364*78e24b28SBarry Smith   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
3656cd56b8bSPeter Brune 
3663542a6bcSPeter Brune   snes->data = (void*) gs;
3673542a6bcSPeter Brune   PetscFunctionReturn(0);
3683542a6bcSPeter Brune }
369