xref: /petsc/src/snes/impls/gs/snesgs.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
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 
23b6266c6eSPeter Brune .seealso: SNESSetTrustRegionTolerance()
24b6266c6eSPeter Brune @*/
25be95d8f1SBarry Smith PetscErrorCode  SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
26b6266c6eSPeter Brune {
27be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
28b6266c6eSPeter Brune 
29b6266c6eSPeter Brune   PetscFunctionBegin;
30b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
31b6266c6eSPeter Brune 
32b6266c6eSPeter Brune   if (abstol != PETSC_DEFAULT) {
33*98921bdaSJacob Faibussowitsch     if (abstol < 0.0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
34b6266c6eSPeter Brune     gs->abstol = abstol;
35b6266c6eSPeter Brune   }
36b6266c6eSPeter Brune   if (rtol != PETSC_DEFAULT) {
37*98921bdaSJacob Faibussowitsch     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ(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   }
40b6266c6eSPeter Brune   if (stol != PETSC_DEFAULT) {
41*98921bdaSJacob Faibussowitsch     if (stol < 0.0) SETERRQ(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) {
45*98921bdaSJacob Faibussowitsch     if (maxit < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
46b6266c6eSPeter Brune     gs->max_its = maxit;
47b6266c6eSPeter Brune   }
48b6266c6eSPeter Brune   PetscFunctionReturn(0);
49b6266c6eSPeter Brune }
50b6266c6eSPeter Brune 
51b6266c6eSPeter Brune /*@
52be95d8f1SBarry Smith    SNESNGSGetTolerances - Gets various parameters used in convergence tests.
53b6266c6eSPeter Brune 
54b6266c6eSPeter Brune    Not Collective
55b6266c6eSPeter Brune 
56b6266c6eSPeter Brune    Input Parameters:
57b6266c6eSPeter Brune +  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 
64b6266c6eSPeter Brune    Notes:
650298fd71SBarry Smith    The user can specify NULL for any parameter that is not needed.
66b6266c6eSPeter Brune 
67b6266c6eSPeter Brune    Level: intermediate
68b6266c6eSPeter Brune 
69b6266c6eSPeter Brune .seealso: SNESSetTolerances()
70b6266c6eSPeter Brune @*/
71be95d8f1SBarry Smith PetscErrorCode  SNESNGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
72b6266c6eSPeter Brune {
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;
81b6266c6eSPeter Brune   PetscFunctionReturn(0);
82b6266c6eSPeter Brune }
836cd56b8bSPeter Brune 
846cd56b8bSPeter Brune /*@
85be95d8f1SBarry Smith    SNESNGSSetSweeps - Sets the number of sweeps of GS to use.
866cd56b8bSPeter Brune 
876cd56b8bSPeter Brune    Input Parameters:
886cd56b8bSPeter Brune +  snes   - the SNES context
896cd56b8bSPeter Brune -  sweeps  - the number of sweeps of GS to perform.
906cd56b8bSPeter Brune 
916cd56b8bSPeter Brune    Level: intermediate
926cd56b8bSPeter Brune 
93d814e60bSBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetNPC(), SNESNGSGetSweeps()
946cd56b8bSPeter Brune @*/
956cd56b8bSPeter Brune 
96be95d8f1SBarry Smith PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
976cd56b8bSPeter Brune {
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 /*@
107be95d8f1SBarry Smith    SNESNGSGetSweeps - Gets the number of sweeps GS will use.
1086cd56b8bSPeter Brune 
1096cd56b8bSPeter Brune    Input Parameters:
1106cd56b8bSPeter Brune .  snes   - the SNES context
1116cd56b8bSPeter Brune 
1126cd56b8bSPeter Brune    Output Parameters:
1136cd56b8bSPeter Brune .  sweeps  - the number of sweeps of GS to perform.
1146cd56b8bSPeter Brune 
1156cd56b8bSPeter Brune    Level: intermediate
1166cd56b8bSPeter Brune 
117d814e60bSBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetNPC(), SNESNGSSetSweeps()
1186cd56b8bSPeter Brune @*/
119be95d8f1SBarry Smith PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt * sweeps)
1206cd56b8bSPeter Brune {
121be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
1226cd56b8bSPeter Brune 
1236cd56b8bSPeter Brune   PetscFunctionBegin;
1246cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1256cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1266cd56b8bSPeter Brune   PetscFunctionReturn(0);
1276cd56b8bSPeter Brune }
1286cd56b8bSPeter Brune 
129be95d8f1SBarry Smith PetscErrorCode SNESReset_NGS(SNES snes)
1303542a6bcSPeter Brune {
1318a86d3c5SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
1328a86d3c5SBarry Smith   PetscErrorCode ierr;
1338a86d3c5SBarry Smith 
1343542a6bcSPeter Brune   PetscFunctionBegin;
1358a86d3c5SBarry Smith   ierr = ISColoringDestroy(&gs->coloring);CHKERRQ(ierr);
1363542a6bcSPeter Brune   PetscFunctionReturn(0);
1373542a6bcSPeter Brune }
1383542a6bcSPeter Brune 
139be95d8f1SBarry Smith PetscErrorCode SNESDestroy_NGS(SNES snes)
1403542a6bcSPeter Brune {
1413542a6bcSPeter Brune   PetscErrorCode ierr;
1423542a6bcSPeter Brune 
1433542a6bcSPeter Brune   PetscFunctionBegin;
144be95d8f1SBarry Smith   ierr = SNESReset_NGS(snes);CHKERRQ(ierr);
14522d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1463542a6bcSPeter Brune   PetscFunctionReturn(0);
1473542a6bcSPeter Brune }
1483542a6bcSPeter Brune 
149be95d8f1SBarry Smith PetscErrorCode SNESSetUp_NGS(SNES snes)
1503542a6bcSPeter Brune {
151ef107cf5SPeter Brune   PetscErrorCode ierr;
152ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
153ef107cf5SPeter Brune 
1543542a6bcSPeter Brune   PetscFunctionBegin;
155be95d8f1SBarry Smith   ierr = SNESGetNGS(snes,&f,NULL);CHKERRQ(ierr);
156ef107cf5SPeter Brune   if (!f) {
157be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
158ef107cf5SPeter Brune   }
1593542a6bcSPeter Brune   PetscFunctionReturn(0);
1603542a6bcSPeter Brune }
1613542a6bcSPeter Brune 
1624416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGS(PetscOptionItems *PetscOptionsObject,SNES snes)
1633542a6bcSPeter Brune {
164be95d8f1SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
1653542a6bcSPeter Brune   PetscErrorCode ierr;
166b6266c6eSPeter Brune   PetscInt       sweeps,max_its=PETSC_DEFAULT;
167b6266c6eSPeter Brune   PetscReal      rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT;
168b6266c6eSPeter Brune   PetscBool      flg,flg1,flg2,flg3;
1693542a6bcSPeter Brune 
1703542a6bcSPeter Brune   PetscFunctionBegin;
171e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES GS options");CHKERRQ(ierr);
1726cd56b8bSPeter Brune   /* GS Options */
173be95d8f1SBarry Smith   ierr = PetscOptionsInt("-snes_ngs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr);
174b6266c6eSPeter Brune   if (flg) {
175be95d8f1SBarry Smith     ierr = SNESNGSSetSweeps(snes,sweeps);CHKERRQ(ierr);
176b6266c6eSPeter Brune   }
17751b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr);
17851b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr);
17951b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr);
18051b6f3e6SMatthew 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);
181b6266c6eSPeter Brune   if (flg || flg1 || flg2 || flg3) {
182be95d8f1SBarry Smith     ierr = SNESNGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr);
183b6266c6eSPeter Brune   }
184ef107cf5SPeter Brune   flg  = PETSC_FALSE;
18598272f38SBarry Smith   ierr = PetscOptionsBool("-snes_ngs_secant","Use finite difference secant approximation with coloring","",flg,&flg,NULL);CHKERRQ(ierr);
186ef107cf5SPeter Brune   if (flg) {
187be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
18898272f38SBarry Smith     ierr = PetscInfo(snes,"Setting default finite difference secant approximation with coloring\n");CHKERRQ(ierr);
189ef107cf5SPeter Brune   }
190be95d8f1SBarry Smith   ierr = PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);CHKERRQ(ierr);
19198272f38SBarry 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);
192ef107cf5SPeter Brune 
1934b32a720SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
1943542a6bcSPeter Brune   PetscFunctionReturn(0);
1953542a6bcSPeter Brune }
1963542a6bcSPeter Brune 
197be95d8f1SBarry Smith PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
1983542a6bcSPeter Brune {
19998272f38SBarry Smith   PetscErrorCode ierr;
20098272f38SBarry Smith   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
20198272f38SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
20298272f38SBarry Smith   PetscBool      iascii;
20398272f38SBarry Smith 
2043542a6bcSPeter Brune   PetscFunctionBegin;
20598272f38SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
20698272f38SBarry Smith   if (iascii) {
20798272f38SBarry Smith     ierr = DMSNESGetNGS(snes->dm,&f,NULL);CHKERRQ(ierr);
20898272f38SBarry Smith     if (f == SNESComputeNGSDefaultSecant) {
209efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Use finite difference secant approximation with coloring with h = %g \n",(double)gs->h);CHKERRQ(ierr);
21098272f38SBarry Smith     }
21198272f38SBarry Smith   }
2123542a6bcSPeter Brune   PetscFunctionReturn(0);
2133542a6bcSPeter Brune }
2143542a6bcSPeter Brune 
215be95d8f1SBarry Smith PetscErrorCode SNESSolve_NGS(SNES snes)
2163542a6bcSPeter Brune {
2173542a6bcSPeter Brune   Vec              F;
2183542a6bcSPeter Brune   Vec              X;
2193542a6bcSPeter Brune   Vec              B;
2203542a6bcSPeter Brune   PetscInt         i;
22106e07b1aSPeter Brune   PetscReal        fnorm;
2223542a6bcSPeter Brune   PetscErrorCode   ierr;
223365a6726SPeter Brune   SNESNormSchedule normschedule;
2243542a6bcSPeter Brune 
2253542a6bcSPeter Brune   PetscFunctionBegin;
226c579b300SPatrick Farrell 
227*98921bdaSJacob Faibussowitsch   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
228c579b300SPatrick Farrell 
229fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
2303542a6bcSPeter Brune   X = snes->vec_sol;
2313542a6bcSPeter Brune   F = snes->vec_func;
2323542a6bcSPeter Brune   B = snes->vec_rhs;
2333542a6bcSPeter Brune 
234e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
235534ebe21SPeter Brune   snes->iter   = 0;
236534ebe21SPeter Brune   snes->norm   = 0.;
237e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
238526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
239534ebe21SPeter Brune 
240365a6726SPeter Brune   ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
241365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2423542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
243e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2443542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2451aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
246e4ed7901SPeter Brune 
2473542a6bcSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
248422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
249e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
250e4ed7901SPeter Brune     snes->iter = 0;
2513542a6bcSPeter Brune     snes->norm = fnorm;
252e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
253a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
254e4ed7901SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
2553542a6bcSPeter Brune 
2563542a6bcSPeter Brune     /* test convergence */
2573542a6bcSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2583542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
259534ebe21SPeter Brune   } else {
260e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
261a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
262534ebe21SPeter Brune   }
263534ebe21SPeter Brune 
2643542a6bcSPeter Brune   /* Call general purpose update function */
2653542a6bcSPeter Brune   if (snes->ops->update) {
2663542a6bcSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2673542a6bcSPeter Brune   }
268534ebe21SPeter Brune 
2693542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
270be95d8f1SBarry Smith     ierr = SNESComputeNGS(snes, B, X);CHKERRQ(ierr);
2714b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
272365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2733542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2743542a6bcSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
275422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
2763542a6bcSPeter Brune       /* Monitor convergence */
277e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
2783542a6bcSPeter Brune       snes->iter = i+1;
2793542a6bcSPeter Brune       snes->norm = fnorm;
280e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
281a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
2823542a6bcSPeter Brune       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
283ac4a6415SMatthew G. Knepley     }
2843542a6bcSPeter Brune     /* Test for convergence */
28560d4fc61SSatish Balay     if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
2863542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
2873542a6bcSPeter Brune     /* Call general purpose update function */
2883542a6bcSPeter Brune     if (snes->ops->update) {
2893542a6bcSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2903542a6bcSPeter Brune     }
2913542a6bcSPeter Brune   }
292365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
293534ebe21SPeter Brune     if (i == snes->max_its) {
294534ebe21SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
295534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
296534ebe21SPeter Brune     }
2971aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
2983542a6bcSPeter Brune   PetscFunctionReturn(0);
2993542a6bcSPeter Brune }
3003542a6bcSPeter Brune 
3013542a6bcSPeter Brune /*MC
30298272f38SBarry Smith   SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation
30398272f38SBarry Smith             using coloring.
3043542a6bcSPeter Brune 
3053542a6bcSPeter Brune    Level: advanced
3063542a6bcSPeter Brune 
30798272f38SBarry Smith   Options Database:
30878e24b28SBarry Smith +   -snes_ngs_sweeps <n> - Number of sweeps of GS to apply
30978e24b28SBarry Smith .    -snes_ngs_atol <atol> - Absolute residual tolerance for GS iteration
31078e24b28SBarry Smith .    -snes_ngs_rtol <rtol> - Relative residual tolerance for GS iteration
31178e24b28SBarry Smith .    -snes_ngs_stol <stol> - Absolute update tolerance for GS iteration
31278e24b28SBarry Smith .    -snes_ngs_max_it <maxit> - Maximum number of sweeps of GS to apply
31398272f38SBarry Smith .    -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is
31498272f38SBarry Smith                         used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring
31598272f38SBarry Smith                         is available or a Jacobian sparse matrix is provided (from which to get the coloring).
31678e24b28SBarry Smith .    -snes_ngs_secant_h <h> - Differencing parameter for secant approximation
31778e24b28SBarry Smith .    -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available.
318a5b23f4aSJose E. Roman -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed
31998272f38SBarry Smith 
3203542a6bcSPeter Brune   Notes:
321d814e60bSBarry Smith   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetNPC(), it will have
3223542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
3233542a6bcSPeter Brune 
32478e24b28SBarry 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
3254f02bc6aSBarry Smith    References:
32696a0c994SBarry Smith .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
3274f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
3284f02bc6aSBarry Smith 
32978e24b28SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetNGS(), SNESType (for list of available types), SNESNGSSetSweeps(), SNESNGSSetTolerances(),
33078e24b28SBarry Smith           SNESSetNormSchedule()
3313542a6bcSPeter Brune M*/
3323542a6bcSPeter Brune 
333be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
3343542a6bcSPeter Brune {
335be95d8f1SBarry Smith   SNES_NGS        *gs;
3363542a6bcSPeter Brune   PetscErrorCode ierr;
3373542a6bcSPeter Brune 
3383542a6bcSPeter Brune   PetscFunctionBegin;
339be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
340be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
341be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
342be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
343be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
344be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3453542a6bcSPeter Brune 
3463542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
347efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
3483542a6bcSPeter Brune 
3494fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_FALSE;
3504fc747eaSLawrence Mitchell 
35188976e71SPeter Brune   if (!snes->tolerancesset) {
3520e444f03SPeter Brune     snes->max_its   = 10000;
3530e444f03SPeter Brune     snes->max_funcs = 10000;
35488976e71SPeter Brune   }
3550e444f03SPeter Brune 
356b00a9115SJed Brown   ierr = PetscNewLog(snes,&gs);CHKERRQ(ierr);
3576cd56b8bSPeter Brune 
3586cd56b8bSPeter Brune   gs->sweeps  = 1;
359b6266c6eSPeter Brune   gs->rtol    = 1e-5;
36078e24b28SBarry Smith   gs->abstol  = PETSC_MACHINE_EPSILON;
36178e24b28SBarry Smith   gs->stol    = 1000*PETSC_MACHINE_EPSILON;
362b6266c6eSPeter Brune   gs->max_its = 50;
36378e24b28SBarry Smith   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
3646cd56b8bSPeter Brune 
3653542a6bcSPeter Brune   snes->data = (void*) gs;
3663542a6bcSPeter Brune   PetscFunctionReturn(0);
3673542a6bcSPeter Brune }
368