xref: /petsc/src/snes/impls/gs/snesgs.c (revision ef107cf556babef7f252fee75ac80f2dee5e62ad)
1b45d2f2cSJed Brown #include <petsc-private/snesimpl.h>             /*I   "petscsnes.h"   I*/
23542a6bcSPeter Brune 
33542a6bcSPeter Brune typedef struct {
4b6266c6eSPeter Brune   PetscInt  sweeps;     /* number of sweeps through the local subdomain before neighbor communication */
5b6266c6eSPeter Brune   PetscInt  max_its;    /* maximum iterations of the inner pointblock solver */
6b6266c6eSPeter Brune   PetscReal rtol;       /* relative tolerance of the inner pointblock solver */
7b6266c6eSPeter Brune   PetscReal abstol;     /* absolute tolerance of the inner pointblock solver */
8b6266c6eSPeter Brune   PetscReal stol;       /* step tolerance of the inner pointblock solver */
93542a6bcSPeter Brune } SNES_GS;
103542a6bcSPeter Brune 
11*ef107cf5SPeter Brune PETSC_EXTERN PetscErrorCode SNESComputeGSDefaultSecant(SNES,Vec,Vec,void *);
126cd56b8bSPeter Brune 
13b6266c6eSPeter Brune #undef __FUNCT__
14b6266c6eSPeter Brune #define __FUNCT__ "SNESGSSetTolerances"
15b6266c6eSPeter Brune /*@
16b6266c6eSPeter Brune    SNESGSSetTolerances - Sets various parameters used in convergence tests.
17b6266c6eSPeter Brune 
18b6266c6eSPeter Brune    Logically Collective on SNES
19b6266c6eSPeter Brune 
20b6266c6eSPeter Brune    Input Parameters:
21b6266c6eSPeter Brune +  snes - the SNES context
22b6266c6eSPeter Brune .  abstol - absolute convergence tolerance
23b6266c6eSPeter Brune .  rtol - relative convergence tolerance
24b6266c6eSPeter Brune .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
25b6266c6eSPeter Brune -  maxit - maximum number of iterations
26b6266c6eSPeter Brune 
27b6266c6eSPeter Brune 
28b6266c6eSPeter Brune    Options Database Keys:
29b6266c6eSPeter Brune +    -snes_gs_atol <abstol> - Sets abstol
30b6266c6eSPeter Brune .    -snes_gs_rtol <rtol> - Sets rtol
31b6266c6eSPeter Brune .    -snes_gs_stol <stol> - Sets stol
32b6266c6eSPeter Brune -    -snes_max_it <maxit> - Sets maxit
33b6266c6eSPeter Brune 
34b6266c6eSPeter Brune    Level: intermediate
35b6266c6eSPeter Brune 
36b6266c6eSPeter Brune .keywords: SNES, nonlinear, gauss-seidel, set, convergence, tolerances
37b6266c6eSPeter Brune 
38b6266c6eSPeter Brune .seealso: SNESSetTrustRegionTolerance()
39b6266c6eSPeter Brune @*/
40b6266c6eSPeter Brune PetscErrorCode  SNESGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
41b6266c6eSPeter Brune {
42b6266c6eSPeter Brune   SNES_GS *gs = (SNES_GS*)snes->data;
43b6266c6eSPeter Brune 
44b6266c6eSPeter Brune   PetscFunctionBegin;
45b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
46b6266c6eSPeter Brune 
47b6266c6eSPeter Brune   if (abstol != PETSC_DEFAULT) {
48ce94432eSBarry Smith     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol);
49b6266c6eSPeter Brune     gs->abstol = abstol;
50b6266c6eSPeter Brune   }
51b6266c6eSPeter Brune   if (rtol != PETSC_DEFAULT) {
52ce94432eSBarry 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",rtol);
53b6266c6eSPeter Brune     gs->rtol = rtol;
54b6266c6eSPeter Brune   }
55b6266c6eSPeter Brune   if (stol != PETSC_DEFAULT) {
56ce94432eSBarry Smith     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %G must be non-negative",stol);
57b6266c6eSPeter Brune     gs->stol = stol;
58b6266c6eSPeter Brune   }
59b6266c6eSPeter Brune   if (maxit != PETSC_DEFAULT) {
60ce94432eSBarry Smith     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
61b6266c6eSPeter Brune     gs->max_its = maxit;
62b6266c6eSPeter Brune   }
63b6266c6eSPeter Brune   PetscFunctionReturn(0);
64b6266c6eSPeter Brune }
65b6266c6eSPeter Brune 
66b6266c6eSPeter Brune #undef __FUNCT__
67b6266c6eSPeter Brune #define __FUNCT__ "SNESGSGetTolerances"
68b6266c6eSPeter Brune /*@
69b6266c6eSPeter Brune    SNESGSGetTolerances - Gets various parameters used in convergence tests.
70b6266c6eSPeter Brune 
71b6266c6eSPeter Brune    Not Collective
72b6266c6eSPeter Brune 
73b6266c6eSPeter Brune    Input Parameters:
74b6266c6eSPeter Brune +  snes - the SNES context
75b6266c6eSPeter Brune .  atol - absolute convergence tolerance
76b6266c6eSPeter Brune .  rtol - relative convergence tolerance
77b6266c6eSPeter Brune .  stol -  convergence tolerance in terms of the norm
78b6266c6eSPeter Brune            of the change in the solution between steps
79b6266c6eSPeter Brune -  maxit - maximum number of iterations
80b6266c6eSPeter Brune 
81b6266c6eSPeter Brune    Notes:
820298fd71SBarry Smith    The user can specify NULL for any parameter that is not needed.
83b6266c6eSPeter Brune 
84b6266c6eSPeter Brune    Level: intermediate
85b6266c6eSPeter Brune 
86b6266c6eSPeter Brune .keywords: SNES, nonlinear, get, convergence, tolerances
87b6266c6eSPeter Brune 
88b6266c6eSPeter Brune .seealso: SNESSetTolerances()
89b6266c6eSPeter Brune @*/
90b6266c6eSPeter Brune PetscErrorCode  SNESGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
91b6266c6eSPeter Brune {
92b6266c6eSPeter Brune   SNES_GS *gs = (SNES_GS*)snes->data;
93b6266c6eSPeter Brune 
94b6266c6eSPeter Brune   PetscFunctionBegin;
95b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
96b6266c6eSPeter Brune   if (atol) *atol = gs->abstol;
97b6266c6eSPeter Brune   if (rtol) *rtol = gs->rtol;
98b6266c6eSPeter Brune   if (stol) *stol = gs->stol;
99b6266c6eSPeter Brune   if (maxit) *maxit = gs->max_its;
100b6266c6eSPeter Brune   PetscFunctionReturn(0);
101b6266c6eSPeter Brune }
1026cd56b8bSPeter Brune 
1036cd56b8bSPeter Brune #undef __FUNCT__
1046cd56b8bSPeter Brune #define __FUNCT__ "SNESGSSetSweeps"
1056cd56b8bSPeter Brune /*@
1066cd56b8bSPeter Brune    SNESGSSetSweeps - Sets the number of sweeps of GS to use.
1076cd56b8bSPeter Brune 
1086cd56b8bSPeter Brune    Input Parameters:
1096cd56b8bSPeter Brune +  snes   - the SNES context
1106cd56b8bSPeter Brune -  sweeps  - the number of sweeps of GS to perform.
1116cd56b8bSPeter Brune 
1126cd56b8bSPeter Brune    Level: intermediate
1136cd56b8bSPeter Brune 
1146cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel
1156cd56b8bSPeter Brune 
1166cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSGetSweeps()
1176cd56b8bSPeter Brune @*/
1186cd56b8bSPeter Brune 
1196cd56b8bSPeter Brune PetscErrorCode SNESGSSetSweeps(SNES snes, PetscInt sweeps)
1206cd56b8bSPeter Brune {
1216cd56b8bSPeter Brune   SNES_GS *gs = (SNES_GS*)snes->data;
1226cd56b8bSPeter Brune 
1236cd56b8bSPeter Brune   PetscFunctionBegin;
1246cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1256cd56b8bSPeter Brune   gs->sweeps = sweeps;
1266cd56b8bSPeter Brune   PetscFunctionReturn(0);
1276cd56b8bSPeter Brune }
1286cd56b8bSPeter Brune 
1296cd56b8bSPeter Brune #undef __FUNCT__
1306cd56b8bSPeter Brune #define __FUNCT__ "SNESGSGetSweeps"
1316cd56b8bSPeter Brune /*@
1326cd56b8bSPeter Brune    SNESGSGetSweeps - Gets the number of sweeps GS will use.
1336cd56b8bSPeter Brune 
1346cd56b8bSPeter Brune    Input Parameters:
1356cd56b8bSPeter Brune .  snes   - the SNES context
1366cd56b8bSPeter Brune 
1376cd56b8bSPeter Brune    Output Parameters:
1386cd56b8bSPeter Brune .  sweeps  - the number of sweeps of GS to perform.
1396cd56b8bSPeter Brune 
1406cd56b8bSPeter Brune    Level: intermediate
1416cd56b8bSPeter Brune 
1426cd56b8bSPeter Brune .keywords: SNES, nonlinear, set, Gauss-Siedel
1436cd56b8bSPeter Brune 
1446cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSSetSweeps()
1456cd56b8bSPeter Brune @*/
1466cd56b8bSPeter Brune PetscErrorCode SNESGSGetSweeps(SNES snes, PetscInt * sweeps)
1476cd56b8bSPeter Brune {
1486cd56b8bSPeter Brune   SNES_GS *gs = (SNES_GS*)snes->data;
1496cd56b8bSPeter Brune 
1506cd56b8bSPeter Brune   PetscFunctionBegin;
1516cd56b8bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1526cd56b8bSPeter Brune   *sweeps = gs->sweeps;
1536cd56b8bSPeter Brune   PetscFunctionReturn(0);
1546cd56b8bSPeter Brune }
1556cd56b8bSPeter Brune 
1566cd56b8bSPeter Brune 
1573542a6bcSPeter Brune #undef __FUNCT__
158c79d6ae8SPeter Brune #define __FUNCT__ "SNESDefaultApplyGS"
1590adebc6cSBarry Smith PetscErrorCode SNESDefaultApplyGS(SNES snes,Vec X,Vec F,void *ctx)
1600adebc6cSBarry Smith {
161c79d6ae8SPeter Brune   PetscFunctionBegin;
162c79d6ae8SPeter Brune   /* see if there's a coloring on the DM */
163c79d6ae8SPeter Brune   PetscFunctionReturn(0);
164c79d6ae8SPeter Brune }
165c79d6ae8SPeter Brune 
166c79d6ae8SPeter Brune #undef __FUNCT__
1673542a6bcSPeter Brune #define __FUNCT__ "SNESReset_GS"
1683542a6bcSPeter Brune PetscErrorCode SNESReset_GS(SNES snes)
1693542a6bcSPeter Brune {
1703542a6bcSPeter Brune   PetscFunctionBegin;
1713542a6bcSPeter Brune   PetscFunctionReturn(0);
1723542a6bcSPeter Brune }
1733542a6bcSPeter Brune 
1743542a6bcSPeter Brune #undef __FUNCT__
1753542a6bcSPeter Brune #define __FUNCT__ "SNESDestroy_GS"
1763542a6bcSPeter Brune PetscErrorCode SNESDestroy_GS(SNES snes)
1773542a6bcSPeter Brune {
1783542a6bcSPeter Brune   PetscErrorCode ierr;
1793542a6bcSPeter Brune 
1803542a6bcSPeter Brune   PetscFunctionBegin;
1813542a6bcSPeter Brune   ierr = SNESReset_GS(snes);CHKERRQ(ierr);
18222d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1833542a6bcSPeter Brune   PetscFunctionReturn(0);
1843542a6bcSPeter Brune }
1853542a6bcSPeter Brune 
1863542a6bcSPeter Brune #undef __FUNCT__
1873542a6bcSPeter Brune #define __FUNCT__ "SNESSetUp_GS"
1883542a6bcSPeter Brune PetscErrorCode SNESSetUp_GS(SNES snes)
1893542a6bcSPeter Brune {
190*ef107cf5SPeter Brune   PetscErrorCode ierr;
191*ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
192*ef107cf5SPeter Brune 
1933542a6bcSPeter Brune   PetscFunctionBegin;
194*ef107cf5SPeter Brune   ierr = SNESGetGS(snes,&f,NULL);CHKERRQ(ierr);
195*ef107cf5SPeter Brune   if (!f) {
196*ef107cf5SPeter Brune     ierr = SNESSetGS(snes,SNESComputeGSDefaultSecant,NULL);CHKERRQ(ierr);
197*ef107cf5SPeter Brune   }
1983542a6bcSPeter Brune   PetscFunctionReturn(0);
1993542a6bcSPeter Brune }
2003542a6bcSPeter Brune 
2013542a6bcSPeter Brune #undef __FUNCT__
2023542a6bcSPeter Brune #define __FUNCT__ "SNESSetFromOptions_GS"
2033542a6bcSPeter Brune PetscErrorCode SNESSetFromOptions_GS(SNES snes)
2043542a6bcSPeter Brune {
2056cd56b8bSPeter Brune   SNES_GS        *gs = (SNES_GS*)snes->data;
2063542a6bcSPeter Brune   PetscErrorCode ierr;
207b6266c6eSPeter Brune   PetscInt       sweeps,max_its=PETSC_DEFAULT;
208b6266c6eSPeter Brune   PetscReal      rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT;
209b6266c6eSPeter Brune   PetscBool      flg,flg1,flg2,flg3;
2103542a6bcSPeter Brune 
2113542a6bcSPeter Brune   PetscFunctionBegin;
2123542a6bcSPeter Brune   ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr);
2136cd56b8bSPeter Brune   /* GS Options */
214b6266c6eSPeter Brune   ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr);
215b6266c6eSPeter Brune   if (flg) {
216b6266c6eSPeter Brune     ierr = SNESGSSetSweeps(snes,sweeps);CHKERRQ(ierr);
217b6266c6eSPeter Brune   }
218b6266c6eSPeter Brune   ierr = PetscOptionsReal("-snes_gs_atol","Number of sweeps of GS to apply","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr);
219b6266c6eSPeter Brune   ierr = PetscOptionsReal("-snes_gs_rtol","Number of sweeps of GS to apply","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr);
220b6266c6eSPeter Brune   ierr = PetscOptionsReal("-snes_gs_stol","Number of sweeps of GS to apply","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr);
221b6266c6eSPeter Brune   ierr = PetscOptionsInt("-snes_gs_max_it","Number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);CHKERRQ(ierr);
222b6266c6eSPeter Brune   if (flg || flg1 || flg2 || flg3) {
223b6266c6eSPeter Brune     ierr = SNESGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr);
224b6266c6eSPeter Brune   }
225*ef107cf5SPeter Brune   flg  = PETSC_FALSE;
226*ef107cf5SPeter Brune   ierr = PetscOptionsBool("-snes_gs_secant","Use pointwise secant local Jacobian approximation","",flg,&flg,NULL);CHKERRQ(ierr);
227*ef107cf5SPeter Brune   if (flg) {
228*ef107cf5SPeter Brune     ierr = SNESSetGS(snes,SNESComputeGSDefaultSecant,NULL);CHKERRQ(ierr);
229*ef107cf5SPeter Brune     ierr = PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");CHKERRQ(ierr);
230*ef107cf5SPeter Brune   }
231*ef107cf5SPeter Brune 
232*ef107cf5SPeter Brune 
2334b32a720SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2343542a6bcSPeter Brune   PetscFunctionReturn(0);
2353542a6bcSPeter Brune }
2363542a6bcSPeter Brune 
2373542a6bcSPeter Brune #undef __FUNCT__
2383542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS"
2393542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer)
2403542a6bcSPeter Brune {
2413542a6bcSPeter Brune   PetscFunctionBegin;
2423542a6bcSPeter Brune   PetscFunctionReturn(0);
2433542a6bcSPeter Brune }
2443542a6bcSPeter Brune 
2453542a6bcSPeter Brune #undef __FUNCT__
2463542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS"
2473542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(SNES snes)
2483542a6bcSPeter Brune {
2493542a6bcSPeter Brune   Vec              F;
2503542a6bcSPeter Brune   Vec              X;
2513542a6bcSPeter Brune   Vec              B;
2523542a6bcSPeter Brune   PetscInt         i;
25306e07b1aSPeter Brune   PetscReal        fnorm;
2543542a6bcSPeter Brune   PetscErrorCode   ierr;
255365a6726SPeter Brune   SNESNormSchedule normschedule;
2563542a6bcSPeter Brune 
2573542a6bcSPeter Brune   PetscFunctionBegin;
2583542a6bcSPeter Brune   X = snes->vec_sol;
2593542a6bcSPeter Brune   F = snes->vec_func;
2603542a6bcSPeter Brune   B = snes->vec_rhs;
2613542a6bcSPeter Brune 
262ce8c27fbSBarry Smith   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
263534ebe21SPeter Brune   snes->iter   = 0;
264534ebe21SPeter Brune   snes->norm   = 0.;
265ce8c27fbSBarry Smith   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
266526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
267534ebe21SPeter Brune 
268365a6726SPeter Brune   ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
269365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2703542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
271e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2723542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2733542a6bcSPeter Brune       if (snes->domainerror) {
2743542a6bcSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2753542a6bcSPeter Brune         PetscFunctionReturn(0);
2763542a6bcSPeter Brune       }
2771aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
278e4ed7901SPeter Brune 
2793542a6bcSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
280189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
281189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
282189a9710SBarry Smith       PetscFunctionReturn(0);
283189a9710SBarry Smith     }
284c1c75074SPeter Brune 
285ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
286e4ed7901SPeter Brune     snes->iter = 0;
2873542a6bcSPeter Brune     snes->norm = fnorm;
288ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
289a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
290e4ed7901SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
2913542a6bcSPeter Brune 
2923542a6bcSPeter Brune     /* test convergence */
2933542a6bcSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2943542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
295534ebe21SPeter Brune   } else {
296ce8c27fbSBarry Smith     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
297a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
298534ebe21SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
299534ebe21SPeter Brune   }
300534ebe21SPeter Brune 
3013542a6bcSPeter Brune   /* Call general purpose update function */
3023542a6bcSPeter Brune   if (snes->ops->update) {
3033542a6bcSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
3043542a6bcSPeter Brune   }
305534ebe21SPeter Brune 
3063542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
3073542a6bcSPeter Brune     ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
3084b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
309365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
3103542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3113542a6bcSPeter Brune       if (snes->domainerror) {
3123542a6bcSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
3133542a6bcSPeter Brune         PetscFunctionReturn(0);
3143542a6bcSPeter Brune       }
3153542a6bcSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
316189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
317189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
318189a9710SBarry Smith         PetscFunctionReturn(0);
319189a9710SBarry Smith       }
320534ebe21SPeter Brune     }
3213542a6bcSPeter Brune     /* Monitor convergence */
322ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3233542a6bcSPeter Brune     snes->iter = i+1;
3243542a6bcSPeter Brune     snes->norm = fnorm;
325ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
326a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
3273542a6bcSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
3283542a6bcSPeter Brune     /* Test for convergence */
329365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3303542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
3313542a6bcSPeter Brune     /* Call general purpose update function */
3323542a6bcSPeter Brune     if (snes->ops->update) {
3333542a6bcSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
3343542a6bcSPeter Brune     }
3353542a6bcSPeter Brune   }
336365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
337534ebe21SPeter Brune     if (i == snes->max_its) {
338534ebe21SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
339534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
340534ebe21SPeter Brune     }
3411aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
3423542a6bcSPeter Brune   PetscFunctionReturn(0);
3433542a6bcSPeter Brune }
3443542a6bcSPeter Brune 
3453542a6bcSPeter Brune /*MC
3463542a6bcSPeter Brune   SNESGS - Just calls the user-provided solution routine provided with SNESSetGS()
3473542a6bcSPeter Brune 
3483542a6bcSPeter Brune    Level: advanced
3493542a6bcSPeter Brune 
3503542a6bcSPeter Brune   Notes:
3513542a6bcSPeter Brune   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetPC(), it will have
3523542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
3533542a6bcSPeter Brune 
3543542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types)
3553542a6bcSPeter Brune M*/
3563542a6bcSPeter Brune 
3573542a6bcSPeter Brune #undef __FUNCT__
3583542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS"
3598cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_GS(SNES snes)
3603542a6bcSPeter Brune {
3613542a6bcSPeter Brune   SNES_GS        *gs;
3623542a6bcSPeter Brune   PetscErrorCode ierr;
3633542a6bcSPeter Brune 
3643542a6bcSPeter Brune   PetscFunctionBegin;
3653542a6bcSPeter Brune   snes->ops->destroy        = SNESDestroy_GS;
3663542a6bcSPeter Brune   snes->ops->setup          = SNESSetUp_GS;
3673542a6bcSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_GS;
3683542a6bcSPeter Brune   snes->ops->view           = SNESView_GS;
3693542a6bcSPeter Brune   snes->ops->solve          = SNESSolve_GS;
3703542a6bcSPeter Brune   snes->ops->reset          = SNESReset_GS;
3713542a6bcSPeter Brune 
3723542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
3733542a6bcSPeter Brune   snes->usespc  = PETSC_FALSE;
3743542a6bcSPeter Brune 
37588976e71SPeter Brune   if (!snes->tolerancesset) {
3760e444f03SPeter Brune     snes->max_its   = 10000;
3770e444f03SPeter Brune     snes->max_funcs = 10000;
37888976e71SPeter Brune   }
3790e444f03SPeter Brune 
3803542a6bcSPeter Brune   ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr);
3816cd56b8bSPeter Brune 
3826cd56b8bSPeter Brune   gs->sweeps  = 1;
383b6266c6eSPeter Brune   gs->rtol    = 1e-5;
384b6266c6eSPeter Brune   gs->abstol  = 1e-15;
385b6266c6eSPeter Brune   gs->stol    = 1e-12;
386b6266c6eSPeter Brune   gs->max_its = 50;
3876cd56b8bSPeter Brune 
3883542a6bcSPeter Brune   snes->data = (void*) gs;
3893542a6bcSPeter Brune   PetscFunctionReturn(0);
3903542a6bcSPeter Brune }
391