xref: /petsc/src/snes/impls/gs/snesgs.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
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 .keywords: SNES, nonlinear, gauss-seidel, set, convergence, tolerances
25b6266c6eSPeter Brune 
26b6266c6eSPeter Brune .seealso: SNESSetTrustRegionTolerance()
27b6266c6eSPeter Brune @*/
28be95d8f1SBarry Smith PetscErrorCode  SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
29b6266c6eSPeter Brune {
30be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
31b6266c6eSPeter Brune 
32b6266c6eSPeter Brune   PetscFunctionBegin;
33b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
34b6266c6eSPeter Brune 
35b6266c6eSPeter Brune   if (abstol != PETSC_DEFAULT) {
3657622a8eSBarry Smith     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
37b6266c6eSPeter Brune     gs->abstol = abstol;
38b6266c6eSPeter Brune   }
39b6266c6eSPeter Brune   if (rtol != PETSC_DEFAULT) {
4057622a8eSBarry 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);
41b6266c6eSPeter Brune     gs->rtol = rtol;
42b6266c6eSPeter Brune   }
43b6266c6eSPeter Brune   if (stol != PETSC_DEFAULT) {
4457622a8eSBarry Smith     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
45b6266c6eSPeter Brune     gs->stol = stol;
46b6266c6eSPeter Brune   }
47b6266c6eSPeter Brune   if (maxit != PETSC_DEFAULT) {
48ce94432eSBarry Smith     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
49b6266c6eSPeter Brune     gs->max_its = maxit;
50b6266c6eSPeter Brune   }
51b6266c6eSPeter Brune   PetscFunctionReturn(0);
52b6266c6eSPeter Brune }
53b6266c6eSPeter Brune 
54b6266c6eSPeter Brune /*@
55be95d8f1SBarry Smith    SNESNGSGetTolerances - Gets various parameters used in convergence tests.
56b6266c6eSPeter Brune 
57b6266c6eSPeter Brune    Not Collective
58b6266c6eSPeter Brune 
59b6266c6eSPeter Brune    Input Parameters:
60b6266c6eSPeter Brune +  snes - the SNES context
61b6266c6eSPeter Brune .  atol - absolute convergence tolerance
62b6266c6eSPeter Brune .  rtol - relative convergence tolerance
63b6266c6eSPeter Brune .  stol -  convergence tolerance in terms of the norm
64b6266c6eSPeter Brune            of the change in the solution between steps
65b6266c6eSPeter Brune -  maxit - maximum number of iterations
66b6266c6eSPeter Brune 
67b6266c6eSPeter Brune    Notes:
680298fd71SBarry Smith    The user can specify NULL for any parameter that is not needed.
69b6266c6eSPeter Brune 
70b6266c6eSPeter Brune    Level: intermediate
71b6266c6eSPeter Brune 
72b6266c6eSPeter Brune .keywords: SNES, nonlinear, get, convergence, tolerances
73b6266c6eSPeter Brune 
74b6266c6eSPeter Brune .seealso: SNESSetTolerances()
75b6266c6eSPeter Brune @*/
76be95d8f1SBarry Smith PetscErrorCode  SNESNGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
77b6266c6eSPeter Brune {
78be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
79b6266c6eSPeter Brune 
80b6266c6eSPeter Brune   PetscFunctionBegin;
81b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
82b6266c6eSPeter Brune   if (atol) *atol = gs->abstol;
83b6266c6eSPeter Brune   if (rtol) *rtol = gs->rtol;
84b6266c6eSPeter Brune   if (stol) *stol = gs->stol;
85b6266c6eSPeter Brune   if (maxit) *maxit = gs->max_its;
86b6266c6eSPeter Brune   PetscFunctionReturn(0);
87b6266c6eSPeter Brune }
886cd56b8bSPeter Brune 
896cd56b8bSPeter Brune /*@
90be95d8f1SBarry Smith    SNESNGSSetSweeps - Sets the number of sweeps of GS to use.
916cd56b8bSPeter Brune 
926cd56b8bSPeter Brune    Input Parameters:
936cd56b8bSPeter Brune +  snes   - the SNES context
946cd56b8bSPeter Brune -  sweeps  - the number of sweeps of GS to perform.
956cd56b8bSPeter Brune 
966cd56b8bSPeter Brune    Level: intermediate
976cd56b8bSPeter Brune 
986cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel
996cd56b8bSPeter Brune 
100be95d8f1SBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetPC(), SNESNGSGetSweeps()
1016cd56b8bSPeter Brune @*/
1026cd56b8bSPeter Brune 
103be95d8f1SBarry Smith PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
1046cd56b8bSPeter Brune {
105be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
1066cd56b8bSPeter Brune 
1076cd56b8bSPeter Brune   PetscFunctionBegin;
1086cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1096cd56b8bSPeter Brune   gs->sweeps = sweeps;
1106cd56b8bSPeter Brune   PetscFunctionReturn(0);
1116cd56b8bSPeter Brune }
1126cd56b8bSPeter Brune 
1136cd56b8bSPeter Brune /*@
114be95d8f1SBarry Smith    SNESNGSGetSweeps - Gets the number of sweeps GS will use.
1156cd56b8bSPeter Brune 
1166cd56b8bSPeter Brune    Input Parameters:
1176cd56b8bSPeter Brune .  snes   - the SNES context
1186cd56b8bSPeter Brune 
1196cd56b8bSPeter Brune    Output Parameters:
1206cd56b8bSPeter Brune .  sweeps  - the number of sweeps of GS to perform.
1216cd56b8bSPeter Brune 
1226cd56b8bSPeter Brune    Level: intermediate
1236cd56b8bSPeter Brune 
1246cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel
1256cd56b8bSPeter Brune 
126be95d8f1SBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetPC(), SNESNGSSetSweeps()
1276cd56b8bSPeter Brune @*/
128be95d8f1SBarry Smith PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt * sweeps)
1296cd56b8bSPeter Brune {
130be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
1316cd56b8bSPeter Brune 
1326cd56b8bSPeter Brune   PetscFunctionBegin;
1336cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1346cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1356cd56b8bSPeter Brune   PetscFunctionReturn(0);
1366cd56b8bSPeter Brune }
1376cd56b8bSPeter Brune 
138be95d8f1SBarry Smith PetscErrorCode SNESReset_NGS(SNES snes)
1393542a6bcSPeter Brune {
1408a86d3c5SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
1418a86d3c5SBarry Smith   PetscErrorCode ierr;
1428a86d3c5SBarry Smith 
1433542a6bcSPeter Brune   PetscFunctionBegin;
1448a86d3c5SBarry Smith   ierr = ISColoringDestroy(&gs->coloring);CHKERRQ(ierr);
1453542a6bcSPeter Brune   PetscFunctionReturn(0);
1463542a6bcSPeter Brune }
1473542a6bcSPeter Brune 
148be95d8f1SBarry Smith PetscErrorCode SNESDestroy_NGS(SNES snes)
1493542a6bcSPeter Brune {
1503542a6bcSPeter Brune   PetscErrorCode ierr;
1513542a6bcSPeter Brune 
1523542a6bcSPeter Brune   PetscFunctionBegin;
153be95d8f1SBarry Smith   ierr = SNESReset_NGS(snes);CHKERRQ(ierr);
15422d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1553542a6bcSPeter Brune   PetscFunctionReturn(0);
1563542a6bcSPeter Brune }
1573542a6bcSPeter Brune 
158be95d8f1SBarry Smith PetscErrorCode SNESSetUp_NGS(SNES snes)
1593542a6bcSPeter Brune {
160ef107cf5SPeter Brune   PetscErrorCode ierr;
161ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
162ef107cf5SPeter Brune 
1633542a6bcSPeter Brune   PetscFunctionBegin;
164be95d8f1SBarry Smith   ierr = SNESGetNGS(snes,&f,NULL);CHKERRQ(ierr);
165ef107cf5SPeter Brune   if (!f) {
166be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
167ef107cf5SPeter Brune   }
1683542a6bcSPeter Brune   PetscFunctionReturn(0);
1693542a6bcSPeter Brune }
1703542a6bcSPeter Brune 
1714416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGS(PetscOptionItems *PetscOptionsObject,SNES snes)
1723542a6bcSPeter Brune {
173be95d8f1SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
1743542a6bcSPeter Brune   PetscErrorCode ierr;
175b6266c6eSPeter Brune   PetscInt       sweeps,max_its=PETSC_DEFAULT;
176b6266c6eSPeter Brune   PetscReal      rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT;
177b6266c6eSPeter Brune   PetscBool      flg,flg1,flg2,flg3;
1783542a6bcSPeter Brune 
1793542a6bcSPeter Brune   PetscFunctionBegin;
180e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES GS options");CHKERRQ(ierr);
1816cd56b8bSPeter Brune   /* GS Options */
182be95d8f1SBarry Smith   ierr = PetscOptionsInt("-snes_ngs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr);
183b6266c6eSPeter Brune   if (flg) {
184be95d8f1SBarry Smith     ierr = SNESNGSSetSweeps(snes,sweeps);CHKERRQ(ierr);
185b6266c6eSPeter Brune   }
18651b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr);
18751b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr);
18851b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr);
18951b6f3e6SMatthew 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);
190b6266c6eSPeter Brune   if (flg || flg1 || flg2 || flg3) {
191be95d8f1SBarry Smith     ierr = SNESNGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr);
192b6266c6eSPeter Brune   }
193ef107cf5SPeter Brune   flg  = PETSC_FALSE;
19498272f38SBarry Smith   ierr = PetscOptionsBool("-snes_ngs_secant","Use finite difference secant approximation with coloring","",flg,&flg,NULL);CHKERRQ(ierr);
195ef107cf5SPeter Brune   if (flg) {
196be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
19798272f38SBarry Smith     ierr = PetscInfo(snes,"Setting default finite difference secant approximation with coloring\n");CHKERRQ(ierr);
198ef107cf5SPeter Brune   }
199be95d8f1SBarry Smith   ierr = PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);CHKERRQ(ierr);
20098272f38SBarry 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);
201ef107cf5SPeter Brune 
2024b32a720SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2033542a6bcSPeter Brune   PetscFunctionReturn(0);
2043542a6bcSPeter Brune }
2053542a6bcSPeter Brune 
206be95d8f1SBarry Smith PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
2073542a6bcSPeter Brune {
20898272f38SBarry Smith   PetscErrorCode ierr;
20998272f38SBarry Smith   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
21098272f38SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
21198272f38SBarry Smith   PetscBool      iascii;
21298272f38SBarry Smith 
2133542a6bcSPeter Brune   PetscFunctionBegin;
21498272f38SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
21598272f38SBarry Smith   if (iascii) {
21698272f38SBarry Smith     ierr = DMSNESGetNGS(snes->dm,&f,NULL);CHKERRQ(ierr);
21798272f38SBarry Smith     if (f == SNESComputeNGSDefaultSecant) {
218*efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Use finite difference secant approximation with coloring with h = %g \n",(double)gs->h);CHKERRQ(ierr);
21998272f38SBarry Smith     }
22098272f38SBarry Smith   }
2213542a6bcSPeter Brune   PetscFunctionReturn(0);
2223542a6bcSPeter Brune }
2233542a6bcSPeter Brune 
224be95d8f1SBarry Smith PetscErrorCode SNESSolve_NGS(SNES snes)
2253542a6bcSPeter Brune {
2263542a6bcSPeter Brune   Vec              F;
2273542a6bcSPeter Brune   Vec              X;
2283542a6bcSPeter Brune   Vec              B;
2293542a6bcSPeter Brune   PetscInt         i;
23006e07b1aSPeter Brune   PetscReal        fnorm;
2313542a6bcSPeter Brune   PetscErrorCode   ierr;
232365a6726SPeter Brune   SNESNormSchedule normschedule;
2333542a6bcSPeter Brune 
2343542a6bcSPeter Brune   PetscFunctionBegin;
235c579b300SPatrick Farrell 
2366c4ed002SBarry 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);
237c579b300SPatrick Farrell 
238fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
2393542a6bcSPeter Brune   X = snes->vec_sol;
2403542a6bcSPeter Brune   F = snes->vec_func;
2413542a6bcSPeter Brune   B = snes->vec_rhs;
2423542a6bcSPeter Brune 
243e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
244534ebe21SPeter Brune   snes->iter   = 0;
245534ebe21SPeter Brune   snes->norm   = 0.;
246e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
247526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
248534ebe21SPeter Brune 
249365a6726SPeter Brune   ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
250365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2513542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
252e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2533542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2541aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
255e4ed7901SPeter Brune 
2563542a6bcSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
257422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
258e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
259e4ed7901SPeter Brune     snes->iter = 0;
2603542a6bcSPeter Brune     snes->norm = fnorm;
261e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
262a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
263e4ed7901SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
2643542a6bcSPeter Brune 
2653542a6bcSPeter Brune     /* test convergence */
2663542a6bcSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2673542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
268534ebe21SPeter Brune   } else {
269e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
270a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
271534ebe21SPeter Brune   }
272534ebe21SPeter Brune 
2733542a6bcSPeter Brune   /* Call general purpose update function */
2743542a6bcSPeter Brune   if (snes->ops->update) {
2753542a6bcSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2763542a6bcSPeter Brune   }
277534ebe21SPeter Brune 
2783542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
279be95d8f1SBarry Smith     ierr = SNESComputeNGS(snes, B, X);CHKERRQ(ierr);
2804b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
281365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2823542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2833542a6bcSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
284422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
2853542a6bcSPeter Brune       /* Monitor convergence */
286e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
2873542a6bcSPeter Brune       snes->iter = i+1;
2883542a6bcSPeter Brune       snes->norm = fnorm;
289e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
290a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
2913542a6bcSPeter Brune       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
292ac4a6415SMatthew G. Knepley     }
2933542a6bcSPeter Brune     /* Test for convergence */
294365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2953542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2963542a6bcSPeter Brune     /* Call general purpose update function */
2973542a6bcSPeter Brune     if (snes->ops->update) {
2983542a6bcSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2993542a6bcSPeter Brune     }
3003542a6bcSPeter Brune   }
301365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
302534ebe21SPeter Brune     if (i == snes->max_its) {
303534ebe21SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
304534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
305534ebe21SPeter Brune     }
3061aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
3073542a6bcSPeter Brune   PetscFunctionReturn(0);
3083542a6bcSPeter Brune }
3093542a6bcSPeter Brune 
3103542a6bcSPeter Brune /*MC
31198272f38SBarry Smith   SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation
31298272f38SBarry Smith             using coloring.
3133542a6bcSPeter Brune 
3143542a6bcSPeter Brune    Level: advanced
3153542a6bcSPeter Brune 
31698272f38SBarry Smith   Options Database:
31798272f38SBarry Smith +   -snes_ngs_sweeps - Number of sweeps of GS to apply
31898272f38SBarry Smith .    -snes_ngs_atol - Absolute residual tolerance for GS iteration
31998272f38SBarry Smith .    -snes_ngs_rtol - Relative residual tolerance for GS iteration
32098272f38SBarry Smith .    -snes_ngs_stol - Absolute update tolerance for GS iteration
32198272f38SBarry Smith .    -snes_ngs_max_it - Maximum number of sweeps of GS to apply
32298272f38SBarry Smith .    -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is
32398272f38SBarry Smith                         used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring
32498272f38SBarry Smith                         is available or a Jacobian sparse matrix is provided (from which to get the coloring).
32598272f38SBarry Smith .    -snes_ngs_secant_h - Differencing parameter for secant approximation
32698272f38SBarry Smith -    -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available.
32798272f38SBarry Smith 
32898272f38SBarry Smith 
3293542a6bcSPeter Brune   Notes:
3303542a6bcSPeter Brune   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetPC(), it will have
3313542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
3323542a6bcSPeter Brune 
3334f02bc6aSBarry Smith    References:
33496a0c994SBarry Smith .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
3354f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
3364f02bc6aSBarry Smith 
33798272f38SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetNGS(), SNESType (for list of available types), SNESNGSSetSweeps(), SNESNGSSetTolerances()
3383542a6bcSPeter Brune M*/
3393542a6bcSPeter Brune 
340be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
3413542a6bcSPeter Brune {
342be95d8f1SBarry Smith   SNES_NGS        *gs;
3433542a6bcSPeter Brune   PetscErrorCode ierr;
3443542a6bcSPeter Brune 
3453542a6bcSPeter Brune   PetscFunctionBegin;
346be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
347be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
348be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
349be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
350be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
351be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3523542a6bcSPeter Brune 
3533542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
354*efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
3553542a6bcSPeter Brune 
3564fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
3574fc747eaSLawrence Mitchell 
35888976e71SPeter Brune   if (!snes->tolerancesset) {
3590e444f03SPeter Brune     snes->max_its   = 10000;
3600e444f03SPeter Brune     snes->max_funcs = 10000;
36188976e71SPeter Brune   }
3620e444f03SPeter Brune 
363b00a9115SJed Brown   ierr = PetscNewLog(snes,&gs);CHKERRQ(ierr);
3646cd56b8bSPeter Brune 
3656cd56b8bSPeter Brune   gs->sweeps  = 1;
366b6266c6eSPeter Brune   gs->rtol    = 1e-5;
367b6266c6eSPeter Brune   gs->abstol  = 1e-15;
368b6266c6eSPeter Brune   gs->stol    = 1e-12;
369b6266c6eSPeter Brune   gs->max_its = 50;
370737a7e12SPeter Brune   gs->h       = 1e-8;
3716cd56b8bSPeter Brune 
3723542a6bcSPeter Brune   snes->data = (void*) gs;
3733542a6bcSPeter Brune   PetscFunctionReturn(0);
3743542a6bcSPeter Brune }
375