xref: /petsc/src/snes/impls/gs/snesgs.c (revision 4f02bc6a7df4e4b03c783003a8e0899ee35fcb05)
1737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h>      /*I "petscsnes.h"  I*/
26cd56b8bSPeter Brune 
3b6266c6eSPeter Brune #undef __FUNCT__
4be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSSetTolerances"
5b6266c6eSPeter Brune /*@
6be95d8f1SBarry Smith    SNESNGSSetTolerances - Sets various parameters used in convergence tests.
7b6266c6eSPeter Brune 
8b6266c6eSPeter Brune    Logically Collective on SNES
9b6266c6eSPeter Brune 
10b6266c6eSPeter Brune    Input Parameters:
11b6266c6eSPeter Brune +  snes - the SNES context
12b6266c6eSPeter Brune .  abstol - absolute convergence tolerance
13b6266c6eSPeter Brune .  rtol - relative convergence tolerance
14b6266c6eSPeter Brune .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
15b6266c6eSPeter Brune -  maxit - maximum number of iterations
16b6266c6eSPeter Brune 
17b6266c6eSPeter Brune 
18b6266c6eSPeter Brune    Options Database Keys:
19be95d8f1SBarry Smith +    -snes_ngs_atol <abstol> - Sets abstol
20be95d8f1SBarry Smith .    -snes_ngs_rtol <rtol> - Sets rtol
21be95d8f1SBarry Smith .    -snes_ngs_stol <stol> - Sets stol
22b6266c6eSPeter Brune -    -snes_max_it <maxit> - Sets maxit
23b6266c6eSPeter Brune 
24b6266c6eSPeter Brune    Level: intermediate
25b6266c6eSPeter Brune 
26b6266c6eSPeter Brune .keywords: SNES, nonlinear, gauss-seidel, set, convergence, tolerances
27b6266c6eSPeter Brune 
28b6266c6eSPeter Brune .seealso: SNESSetTrustRegionTolerance()
29b6266c6eSPeter Brune @*/
30be95d8f1SBarry Smith PetscErrorCode  SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
31b6266c6eSPeter Brune {
32be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
33b6266c6eSPeter Brune 
34b6266c6eSPeter Brune   PetscFunctionBegin;
35b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
36b6266c6eSPeter Brune 
37b6266c6eSPeter Brune   if (abstol != PETSC_DEFAULT) {
3857622a8eSBarry Smith     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
39b6266c6eSPeter Brune     gs->abstol = abstol;
40b6266c6eSPeter Brune   }
41b6266c6eSPeter Brune   if (rtol != PETSC_DEFAULT) {
4257622a8eSBarry 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);
43b6266c6eSPeter Brune     gs->rtol = rtol;
44b6266c6eSPeter Brune   }
45b6266c6eSPeter Brune   if (stol != PETSC_DEFAULT) {
4657622a8eSBarry Smith     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
47b6266c6eSPeter Brune     gs->stol = stol;
48b6266c6eSPeter Brune   }
49b6266c6eSPeter Brune   if (maxit != PETSC_DEFAULT) {
50ce94432eSBarry Smith     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
51b6266c6eSPeter Brune     gs->max_its = maxit;
52b6266c6eSPeter Brune   }
53b6266c6eSPeter Brune   PetscFunctionReturn(0);
54b6266c6eSPeter Brune }
55b6266c6eSPeter Brune 
56b6266c6eSPeter Brune #undef __FUNCT__
57be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSGetTolerances"
58b6266c6eSPeter Brune /*@
59be95d8f1SBarry Smith    SNESNGSGetTolerances - Gets various parameters used in convergence tests.
60b6266c6eSPeter Brune 
61b6266c6eSPeter Brune    Not Collective
62b6266c6eSPeter Brune 
63b6266c6eSPeter Brune    Input Parameters:
64b6266c6eSPeter Brune +  snes - the SNES context
65b6266c6eSPeter Brune .  atol - absolute convergence tolerance
66b6266c6eSPeter Brune .  rtol - relative convergence tolerance
67b6266c6eSPeter Brune .  stol -  convergence tolerance in terms of the norm
68b6266c6eSPeter Brune            of the change in the solution between steps
69b6266c6eSPeter Brune -  maxit - maximum number of iterations
70b6266c6eSPeter Brune 
71b6266c6eSPeter Brune    Notes:
720298fd71SBarry Smith    The user can specify NULL for any parameter that is not needed.
73b6266c6eSPeter Brune 
74b6266c6eSPeter Brune    Level: intermediate
75b6266c6eSPeter Brune 
76b6266c6eSPeter Brune .keywords: SNES, nonlinear, get, convergence, tolerances
77b6266c6eSPeter Brune 
78b6266c6eSPeter Brune .seealso: SNESSetTolerances()
79b6266c6eSPeter Brune @*/
80be95d8f1SBarry Smith PetscErrorCode  SNESNGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
81b6266c6eSPeter Brune {
82be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
83b6266c6eSPeter Brune 
84b6266c6eSPeter Brune   PetscFunctionBegin;
85b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
86b6266c6eSPeter Brune   if (atol) *atol = gs->abstol;
87b6266c6eSPeter Brune   if (rtol) *rtol = gs->rtol;
88b6266c6eSPeter Brune   if (stol) *stol = gs->stol;
89b6266c6eSPeter Brune   if (maxit) *maxit = gs->max_its;
90b6266c6eSPeter Brune   PetscFunctionReturn(0);
91b6266c6eSPeter Brune }
926cd56b8bSPeter Brune 
936cd56b8bSPeter Brune #undef __FUNCT__
94be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSSetSweeps"
956cd56b8bSPeter Brune /*@
96be95d8f1SBarry Smith    SNESNGSSetSweeps - Sets the number of sweeps of GS to use.
976cd56b8bSPeter Brune 
986cd56b8bSPeter Brune    Input Parameters:
996cd56b8bSPeter Brune +  snes   - the SNES context
1006cd56b8bSPeter Brune -  sweeps  - the number of sweeps of GS to perform.
1016cd56b8bSPeter Brune 
1026cd56b8bSPeter Brune    Level: intermediate
1036cd56b8bSPeter Brune 
1046cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel
1056cd56b8bSPeter Brune 
106be95d8f1SBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetPC(), SNESNGSGetSweeps()
1076cd56b8bSPeter Brune @*/
1086cd56b8bSPeter Brune 
109be95d8f1SBarry Smith PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
1106cd56b8bSPeter Brune {
111be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
1126cd56b8bSPeter Brune 
1136cd56b8bSPeter Brune   PetscFunctionBegin;
1146cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1156cd56b8bSPeter Brune   gs->sweeps = sweeps;
1166cd56b8bSPeter Brune   PetscFunctionReturn(0);
1176cd56b8bSPeter Brune }
1186cd56b8bSPeter Brune 
1196cd56b8bSPeter Brune #undef __FUNCT__
120be95d8f1SBarry Smith #define __FUNCT__ "SNESNGSGetSweeps"
1216cd56b8bSPeter Brune /*@
122be95d8f1SBarry Smith    SNESNGSGetSweeps - Gets the number of sweeps GS will use.
1236cd56b8bSPeter Brune 
1246cd56b8bSPeter Brune    Input Parameters:
1256cd56b8bSPeter Brune .  snes   - the SNES context
1266cd56b8bSPeter Brune 
1276cd56b8bSPeter Brune    Output Parameters:
1286cd56b8bSPeter Brune .  sweeps  - the number of sweeps of GS to perform.
1296cd56b8bSPeter Brune 
1306cd56b8bSPeter Brune    Level: intermediate
1316cd56b8bSPeter Brune 
1326cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel
1336cd56b8bSPeter Brune 
134be95d8f1SBarry Smith .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetPC(), SNESNGSSetSweeps()
1356cd56b8bSPeter Brune @*/
136be95d8f1SBarry Smith PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt * sweeps)
1376cd56b8bSPeter Brune {
138be95d8f1SBarry Smith   SNES_NGS *gs = (SNES_NGS*)snes->data;
1396cd56b8bSPeter Brune 
1406cd56b8bSPeter Brune   PetscFunctionBegin;
1416cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1426cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1436cd56b8bSPeter Brune   PetscFunctionReturn(0);
1446cd56b8bSPeter Brune }
1456cd56b8bSPeter Brune 
1466cd56b8bSPeter Brune 
1473542a6bcSPeter Brune #undef __FUNCT__
148be95d8f1SBarry Smith #define __FUNCT__ "SNESDefaultApplyNGS"
149be95d8f1SBarry Smith PetscErrorCode SNESDefaultApplyNGS(SNES snes,Vec X,Vec F,void *ctx)
1500adebc6cSBarry Smith {
151c79d6ae8SPeter Brune   PetscFunctionBegin;
152c79d6ae8SPeter Brune   /* see if there's a coloring on the DM */
153c79d6ae8SPeter Brune   PetscFunctionReturn(0);
154c79d6ae8SPeter Brune }
155c79d6ae8SPeter Brune 
156c79d6ae8SPeter Brune #undef __FUNCT__
157be95d8f1SBarry Smith #define __FUNCT__ "SNESReset_NGS"
158be95d8f1SBarry Smith PetscErrorCode SNESReset_NGS(SNES snes)
1593542a6bcSPeter Brune {
1603542a6bcSPeter Brune   PetscFunctionBegin;
1613542a6bcSPeter Brune   PetscFunctionReturn(0);
1623542a6bcSPeter Brune }
1633542a6bcSPeter Brune 
1643542a6bcSPeter Brune #undef __FUNCT__
165be95d8f1SBarry Smith #define __FUNCT__ "SNESDestroy_NGS"
166be95d8f1SBarry Smith PetscErrorCode SNESDestroy_NGS(SNES snes)
1673542a6bcSPeter Brune {
1683542a6bcSPeter Brune   PetscErrorCode ierr;
1693542a6bcSPeter Brune 
1703542a6bcSPeter Brune   PetscFunctionBegin;
171be95d8f1SBarry Smith   ierr = SNESReset_NGS(snes);CHKERRQ(ierr);
17222d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1733542a6bcSPeter Brune   PetscFunctionReturn(0);
1743542a6bcSPeter Brune }
1753542a6bcSPeter Brune 
1763542a6bcSPeter Brune #undef __FUNCT__
177be95d8f1SBarry Smith #define __FUNCT__ "SNESSetUp_NGS"
178be95d8f1SBarry Smith PetscErrorCode SNESSetUp_NGS(SNES snes)
1793542a6bcSPeter Brune {
180ef107cf5SPeter Brune   PetscErrorCode ierr;
181ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
182ef107cf5SPeter Brune 
1833542a6bcSPeter Brune   PetscFunctionBegin;
184be95d8f1SBarry Smith   ierr = SNESGetNGS(snes,&f,NULL);CHKERRQ(ierr);
185ef107cf5SPeter Brune   if (!f) {
186be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
187ef107cf5SPeter Brune   }
1883542a6bcSPeter Brune   PetscFunctionReturn(0);
1893542a6bcSPeter Brune }
1903542a6bcSPeter Brune 
1913542a6bcSPeter Brune #undef __FUNCT__
192be95d8f1SBarry Smith #define __FUNCT__ "SNESSetFromOptions_NGS"
1934416b707SBarry Smith PetscErrorCode SNESSetFromOptions_NGS(PetscOptionItems *PetscOptionsObject,SNES snes)
1943542a6bcSPeter Brune {
195be95d8f1SBarry Smith   SNES_NGS       *gs = (SNES_NGS*)snes->data;
1963542a6bcSPeter Brune   PetscErrorCode ierr;
197b6266c6eSPeter Brune   PetscInt       sweeps,max_its=PETSC_DEFAULT;
198b6266c6eSPeter Brune   PetscReal      rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT;
199b6266c6eSPeter Brune   PetscBool      flg,flg1,flg2,flg3;
2003542a6bcSPeter Brune 
2013542a6bcSPeter Brune   PetscFunctionBegin;
202e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES GS options");CHKERRQ(ierr);
2036cd56b8bSPeter Brune   /* GS Options */
204be95d8f1SBarry Smith   ierr = PetscOptionsInt("-snes_ngs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr);
205b6266c6eSPeter Brune   if (flg) {
206be95d8f1SBarry Smith     ierr = SNESNGSSetSweeps(snes,sweeps);CHKERRQ(ierr);
207b6266c6eSPeter Brune   }
20851b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr);
20951b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr);
21051b6f3e6SMatthew G. Knepley   ierr = PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr);
21151b6f3e6SMatthew 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);
212b6266c6eSPeter Brune   if (flg || flg1 || flg2 || flg3) {
213be95d8f1SBarry Smith     ierr = SNESNGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr);
214b6266c6eSPeter Brune   }
215ef107cf5SPeter Brune   flg  = PETSC_FALSE;
216be95d8f1SBarry Smith   ierr = PetscOptionsBool("-snes_ngs_secant","Use pointwise secant local Jacobian approximation","",flg,&flg,NULL);CHKERRQ(ierr);
217ef107cf5SPeter Brune   if (flg) {
218be95d8f1SBarry Smith     ierr = SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);CHKERRQ(ierr);
219ef107cf5SPeter Brune     ierr = PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");CHKERRQ(ierr);
220ef107cf5SPeter Brune   }
221be95d8f1SBarry Smith   ierr = PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);CHKERRQ(ierr);
222be95d8f1SBarry Smith   ierr = PetscOptionsBool("-snes_ngs_secant_mat_coloring","Use the Jacobian coloring for the secant GS","",gs->secant_mat,&gs->secant_mat,&flg);CHKERRQ(ierr);
223ef107cf5SPeter Brune 
2244b32a720SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2253542a6bcSPeter Brune   PetscFunctionReturn(0);
2263542a6bcSPeter Brune }
2273542a6bcSPeter Brune 
2283542a6bcSPeter Brune #undef __FUNCT__
229be95d8f1SBarry Smith #define __FUNCT__ "SNESView_NGS"
230be95d8f1SBarry Smith PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
2313542a6bcSPeter Brune {
2323542a6bcSPeter Brune   PetscFunctionBegin;
2333542a6bcSPeter Brune   PetscFunctionReturn(0);
2343542a6bcSPeter Brune }
2353542a6bcSPeter Brune 
2363542a6bcSPeter Brune #undef __FUNCT__
237be95d8f1SBarry Smith #define __FUNCT__ "SNESSolve_NGS"
238be95d8f1SBarry Smith PetscErrorCode SNESSolve_NGS(SNES snes)
2393542a6bcSPeter Brune {
2403542a6bcSPeter Brune   Vec              F;
2413542a6bcSPeter Brune   Vec              X;
2423542a6bcSPeter Brune   Vec              B;
2433542a6bcSPeter Brune   PetscInt         i;
24406e07b1aSPeter Brune   PetscReal        fnorm;
2453542a6bcSPeter Brune   PetscErrorCode   ierr;
246365a6726SPeter Brune   SNESNormSchedule normschedule;
2473542a6bcSPeter Brune 
2483542a6bcSPeter Brune   PetscFunctionBegin;
249c579b300SPatrick Farrell 
2506c4ed002SBarry 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);
251c579b300SPatrick Farrell 
252fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
2533542a6bcSPeter Brune   X = snes->vec_sol;
2543542a6bcSPeter Brune   F = snes->vec_func;
2553542a6bcSPeter Brune   B = snes->vec_rhs;
2563542a6bcSPeter Brune 
257e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
258534ebe21SPeter Brune   snes->iter   = 0;
259534ebe21SPeter Brune   snes->norm   = 0.;
260e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
261526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
262534ebe21SPeter Brune 
263365a6726SPeter Brune   ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
264365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2653542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
266e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2673542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2681aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
269e4ed7901SPeter Brune 
2703542a6bcSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
271422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
272e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
273e4ed7901SPeter Brune     snes->iter = 0;
2743542a6bcSPeter Brune     snes->norm = fnorm;
275e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
276a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
277e4ed7901SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
2783542a6bcSPeter Brune 
2793542a6bcSPeter Brune     /* test convergence */
2803542a6bcSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2813542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
282534ebe21SPeter Brune   } else {
283e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
284a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
285534ebe21SPeter Brune   }
286534ebe21SPeter Brune 
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   }
291534ebe21SPeter Brune 
2923542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
293be95d8f1SBarry Smith     ierr = SNESComputeNGS(snes, B, X);CHKERRQ(ierr);
2944b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
295365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
2963542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2973542a6bcSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
298422a814eSBarry Smith       SNESCheckFunctionNorm(snes,fnorm);
2993542a6bcSPeter Brune       /* Monitor convergence */
300e04113cfSBarry Smith       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3013542a6bcSPeter Brune       snes->iter = i+1;
3023542a6bcSPeter Brune       snes->norm = fnorm;
303e04113cfSBarry Smith       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
304a71f0d7dSBarry Smith       ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
3053542a6bcSPeter Brune       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
306ac4a6415SMatthew G. Knepley     }
3073542a6bcSPeter Brune     /* Test for convergence */
308365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3093542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
3103542a6bcSPeter Brune     /* Call general purpose update function */
3113542a6bcSPeter Brune     if (snes->ops->update) {
3123542a6bcSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
3133542a6bcSPeter Brune     }
3143542a6bcSPeter Brune   }
315365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
316534ebe21SPeter Brune     if (i == snes->max_its) {
317534ebe21SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
318534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
319534ebe21SPeter Brune     }
3201aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
3213542a6bcSPeter Brune   PetscFunctionReturn(0);
3223542a6bcSPeter Brune }
3233542a6bcSPeter Brune 
3243542a6bcSPeter Brune /*MC
325be95d8f1SBarry Smith   SNESNGS - Just calls the user-provided solution routine provided with SNESSetNGS()
3263542a6bcSPeter Brune 
3273542a6bcSPeter Brune    Level: advanced
3283542a6bcSPeter Brune 
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 
333*4f02bc6aSBarry Smith    References:
334*4f02bc6aSBarry Smith 
335*4f02bc6aSBarry Smith    Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
336*4f02bc6aSBarry Smith    SIAM Review, 57(4), 2015
337*4f02bc6aSBarry Smith 
338be95d8f1SBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetNGS(), SNESType (for list of available types)
3393542a6bcSPeter Brune M*/
3403542a6bcSPeter Brune 
3413542a6bcSPeter Brune #undef __FUNCT__
342be95d8f1SBarry Smith #define __FUNCT__ "SNESCreate_NGS"
343be95d8f1SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
3443542a6bcSPeter Brune {
345be95d8f1SBarry Smith   SNES_NGS        *gs;
3463542a6bcSPeter Brune   PetscErrorCode ierr;
3473542a6bcSPeter Brune 
3483542a6bcSPeter Brune   PetscFunctionBegin;
349be95d8f1SBarry Smith   snes->ops->destroy        = SNESDestroy_NGS;
350be95d8f1SBarry Smith   snes->ops->setup          = SNESSetUp_NGS;
351be95d8f1SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
352be95d8f1SBarry Smith   snes->ops->view           = SNESView_NGS;
353be95d8f1SBarry Smith   snes->ops->solve          = SNESSolve_NGS;
354be95d8f1SBarry Smith   snes->ops->reset          = SNESReset_NGS;
3553542a6bcSPeter Brune 
3563542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
3573542a6bcSPeter Brune   snes->usespc  = PETSC_FALSE;
3583542a6bcSPeter Brune 
35988976e71SPeter Brune   if (!snes->tolerancesset) {
3600e444f03SPeter Brune     snes->max_its   = 10000;
3610e444f03SPeter Brune     snes->max_funcs = 10000;
36288976e71SPeter Brune   }
3630e444f03SPeter Brune 
364b00a9115SJed Brown   ierr = PetscNewLog(snes,&gs);CHKERRQ(ierr);
3656cd56b8bSPeter Brune 
3666cd56b8bSPeter Brune   gs->sweeps  = 1;
367b6266c6eSPeter Brune   gs->rtol    = 1e-5;
368b6266c6eSPeter Brune   gs->abstol  = 1e-15;
369b6266c6eSPeter Brune   gs->stol    = 1e-12;
370b6266c6eSPeter Brune   gs->max_its = 50;
371737a7e12SPeter Brune   gs->h       = 1e-8;
3726cd56b8bSPeter Brune 
3733542a6bcSPeter Brune   snes->data = (void*) gs;
3743542a6bcSPeter Brune   PetscFunctionReturn(0);
3753542a6bcSPeter Brune }
376