xref: /petsc/src/snes/impls/gs/snesgs.c (revision 737a7e128d3ed57ea128fd0d033fca43bac7efb0)
1*737a7e12SPeter Brune #include <../src/snes/impls/gs/gsimpl.h>      /*I "petscsnes.h"  I*/
26cd56b8bSPeter Brune 
3b6266c6eSPeter Brune #undef __FUNCT__
4b6266c6eSPeter Brune #define __FUNCT__ "SNESGSSetTolerances"
5b6266c6eSPeter Brune /*@
6b6266c6eSPeter Brune    SNESGSSetTolerances - 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:
19b6266c6eSPeter Brune +    -snes_gs_atol <abstol> - Sets abstol
20b6266c6eSPeter Brune .    -snes_gs_rtol <rtol> - Sets rtol
21b6266c6eSPeter Brune .    -snes_gs_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 @*/
30b6266c6eSPeter Brune PetscErrorCode  SNESGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
31b6266c6eSPeter Brune {
32b6266c6eSPeter Brune   SNES_GS *gs = (SNES_GS*)snes->data;
33b6266c6eSPeter Brune 
34b6266c6eSPeter Brune   PetscFunctionBegin;
35b6266c6eSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
36b6266c6eSPeter Brune 
37b6266c6eSPeter Brune   if (abstol != PETSC_DEFAULT) {
38ce94432eSBarry Smith     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol);
39b6266c6eSPeter Brune     gs->abstol = abstol;
40b6266c6eSPeter Brune   }
41b6266c6eSPeter Brune   if (rtol != PETSC_DEFAULT) {
42ce94432eSBarry 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);
43b6266c6eSPeter Brune     gs->rtol = rtol;
44b6266c6eSPeter Brune   }
45b6266c6eSPeter Brune   if (stol != PETSC_DEFAULT) {
46ce94432eSBarry Smith     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %G must be non-negative",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__
57b6266c6eSPeter Brune #define __FUNCT__ "SNESGSGetTolerances"
58b6266c6eSPeter Brune /*@
59b6266c6eSPeter Brune    SNESGSGetTolerances - 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 @*/
80b6266c6eSPeter Brune PetscErrorCode  SNESGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
81b6266c6eSPeter Brune {
82b6266c6eSPeter Brune   SNES_GS *gs = (SNES_GS*)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__
946cd56b8bSPeter Brune #define __FUNCT__ "SNESGSSetSweeps"
956cd56b8bSPeter Brune /*@
966cd56b8bSPeter Brune    SNESGSSetSweeps - 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 
1066cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSGetSweeps()
1076cd56b8bSPeter Brune @*/
1086cd56b8bSPeter Brune 
1096cd56b8bSPeter Brune PetscErrorCode SNESGSSetSweeps(SNES snes, PetscInt sweeps)
1106cd56b8bSPeter Brune {
1116cd56b8bSPeter Brune   SNES_GS *gs = (SNES_GS*)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__
1206cd56b8bSPeter Brune #define __FUNCT__ "SNESGSGetSweeps"
1216cd56b8bSPeter Brune /*@
1226cd56b8bSPeter Brune    SNESGSGetSweeps - 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 
1346cd56b8bSPeter Brune .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSSetSweeps()
1356cd56b8bSPeter Brune @*/
1366cd56b8bSPeter Brune PetscErrorCode SNESGSGetSweeps(SNES snes, PetscInt * sweeps)
1376cd56b8bSPeter Brune {
1386cd56b8bSPeter Brune   SNES_GS *gs = (SNES_GS*)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__
148c79d6ae8SPeter Brune #define __FUNCT__ "SNESDefaultApplyGS"
1490adebc6cSBarry Smith PetscErrorCode SNESDefaultApplyGS(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__
1573542a6bcSPeter Brune #define __FUNCT__ "SNESReset_GS"
1583542a6bcSPeter Brune PetscErrorCode SNESReset_GS(SNES snes)
1593542a6bcSPeter Brune {
1603542a6bcSPeter Brune   PetscFunctionBegin;
1613542a6bcSPeter Brune   PetscFunctionReturn(0);
1623542a6bcSPeter Brune }
1633542a6bcSPeter Brune 
1643542a6bcSPeter Brune #undef __FUNCT__
1653542a6bcSPeter Brune #define __FUNCT__ "SNESDestroy_GS"
1663542a6bcSPeter Brune PetscErrorCode SNESDestroy_GS(SNES snes)
1673542a6bcSPeter Brune {
1683542a6bcSPeter Brune   PetscErrorCode ierr;
1693542a6bcSPeter Brune 
1703542a6bcSPeter Brune   PetscFunctionBegin;
1713542a6bcSPeter Brune   ierr = SNESReset_GS(snes);CHKERRQ(ierr);
17222d28d08SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
1733542a6bcSPeter Brune   PetscFunctionReturn(0);
1743542a6bcSPeter Brune }
1753542a6bcSPeter Brune 
1763542a6bcSPeter Brune #undef __FUNCT__
1773542a6bcSPeter Brune #define __FUNCT__ "SNESSetUp_GS"
1783542a6bcSPeter Brune PetscErrorCode SNESSetUp_GS(SNES snes)
1793542a6bcSPeter Brune {
180ef107cf5SPeter Brune   PetscErrorCode ierr;
181ef107cf5SPeter Brune   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
182ef107cf5SPeter Brune 
1833542a6bcSPeter Brune   PetscFunctionBegin;
184ef107cf5SPeter Brune   ierr = SNESGetGS(snes,&f,NULL);CHKERRQ(ierr);
185ef107cf5SPeter Brune   if (!f) {
186ef107cf5SPeter Brune     ierr = SNESSetGS(snes,SNESComputeGSDefaultSecant,NULL);CHKERRQ(ierr);
187ef107cf5SPeter Brune   }
1883542a6bcSPeter Brune   PetscFunctionReturn(0);
1893542a6bcSPeter Brune }
1903542a6bcSPeter Brune 
1913542a6bcSPeter Brune #undef __FUNCT__
1923542a6bcSPeter Brune #define __FUNCT__ "SNESSetFromOptions_GS"
1933542a6bcSPeter Brune PetscErrorCode SNESSetFromOptions_GS(SNES snes)
1943542a6bcSPeter Brune {
1956cd56b8bSPeter Brune   SNES_GS        *gs = (SNES_GS*)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;
2023542a6bcSPeter Brune   ierr = PetscOptionsHead("SNES GS options");CHKERRQ(ierr);
2036cd56b8bSPeter Brune   /* GS Options */
204b6266c6eSPeter Brune   ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);CHKERRQ(ierr);
205b6266c6eSPeter Brune   if (flg) {
206b6266c6eSPeter Brune     ierr = SNESGSSetSweeps(snes,sweeps);CHKERRQ(ierr);
207b6266c6eSPeter Brune   }
208b6266c6eSPeter Brune   ierr = PetscOptionsReal("-snes_gs_atol","Number of sweeps of GS to apply","SNESComputeGS",gs->abstol,&atol,&flg);CHKERRQ(ierr);
209b6266c6eSPeter Brune   ierr = PetscOptionsReal("-snes_gs_rtol","Number of sweeps of GS to apply","SNESComputeGS",gs->rtol,&rtol,&flg1);CHKERRQ(ierr);
210b6266c6eSPeter Brune   ierr = PetscOptionsReal("-snes_gs_stol","Number of sweeps of GS to apply","SNESComputeGS",gs->stol,&stol,&flg2);CHKERRQ(ierr);
211b6266c6eSPeter Brune   ierr = PetscOptionsInt("-snes_gs_max_it","Number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);CHKERRQ(ierr);
212b6266c6eSPeter Brune   if (flg || flg1 || flg2 || flg3) {
213b6266c6eSPeter Brune     ierr = SNESGSSetTolerances(snes,atol,rtol,stol,max_its);CHKERRQ(ierr);
214b6266c6eSPeter Brune   }
215ef107cf5SPeter Brune   flg  = PETSC_FALSE;
216ef107cf5SPeter Brune   ierr = PetscOptionsBool("-snes_gs_secant","Use pointwise secant local Jacobian approximation","",flg,&flg,NULL);CHKERRQ(ierr);
217ef107cf5SPeter Brune   if (flg) {
218ef107cf5SPeter Brune     ierr = SNESSetGS(snes,SNESComputeGSDefaultSecant,NULL);CHKERRQ(ierr);
219ef107cf5SPeter Brune     ierr = PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");CHKERRQ(ierr);
220ef107cf5SPeter Brune   }
221*737a7e12SPeter Brune   ierr = PetscOptionsReal("-snes_gs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);CHKERRQ(ierr);
222*737a7e12SPeter Brune   ierr = PetscOptionsBool("-snes_gs_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__
2293542a6bcSPeter Brune #define __FUNCT__ "SNESView_GS"
2303542a6bcSPeter Brune PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer)
2313542a6bcSPeter Brune {
2323542a6bcSPeter Brune   PetscFunctionBegin;
2333542a6bcSPeter Brune   PetscFunctionReturn(0);
2343542a6bcSPeter Brune }
2353542a6bcSPeter Brune 
2363542a6bcSPeter Brune #undef __FUNCT__
2373542a6bcSPeter Brune #define __FUNCT__ "SNESSolve_GS"
2383542a6bcSPeter Brune PetscErrorCode SNESSolve_GS(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;
2493542a6bcSPeter Brune   X = snes->vec_sol;
2503542a6bcSPeter Brune   F = snes->vec_func;
2513542a6bcSPeter Brune   B = snes->vec_rhs;
2523542a6bcSPeter Brune 
253ce8c27fbSBarry Smith   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
254534ebe21SPeter Brune   snes->iter   = 0;
255534ebe21SPeter Brune   snes->norm   = 0.;
256ce8c27fbSBarry Smith   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
257526b802eSJed Brown   snes->reason = SNES_CONVERGED_ITERATING;
258534ebe21SPeter Brune 
259365a6726SPeter Brune   ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
260365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
2613542a6bcSPeter Brune     /* compute the initial function and preconditioned update delX */
262e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
2633542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2643542a6bcSPeter Brune       if (snes->domainerror) {
2653542a6bcSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2663542a6bcSPeter Brune         PetscFunctionReturn(0);
2673542a6bcSPeter Brune       }
2681aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
269e4ed7901SPeter Brune 
2703542a6bcSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
271189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
272189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
273189a9710SBarry Smith       PetscFunctionReturn(0);
274189a9710SBarry Smith     }
275c1c75074SPeter Brune 
276ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
277e4ed7901SPeter Brune     snes->iter = 0;
2783542a6bcSPeter Brune     snes->norm = fnorm;
279ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
280a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
281e4ed7901SPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
2823542a6bcSPeter Brune 
2833542a6bcSPeter Brune     /* test convergence */
2843542a6bcSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2853542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
286534ebe21SPeter Brune   } else {
287ce8c27fbSBarry Smith     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
288a71f0d7dSBarry Smith     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
289534ebe21SPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
290534ebe21SPeter Brune   }
291534ebe21SPeter Brune 
2923542a6bcSPeter Brune   /* Call general purpose update function */
2933542a6bcSPeter Brune   if (snes->ops->update) {
2943542a6bcSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
2953542a6bcSPeter Brune   }
296534ebe21SPeter Brune 
2973542a6bcSPeter Brune   for (i = 0; i < snes->max_its; i++) {
2983542a6bcSPeter Brune     ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
2994b32a720SPeter Brune     /* only compute norms if requested or about to exit due to maximum iterations */
300365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
3013542a6bcSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
3023542a6bcSPeter Brune       if (snes->domainerror) {
3033542a6bcSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
3043542a6bcSPeter Brune         PetscFunctionReturn(0);
3053542a6bcSPeter Brune       }
3063542a6bcSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
307189a9710SBarry Smith       if (PetscIsInfOrNanReal(fnorm)) {
308189a9710SBarry Smith         snes->reason = SNES_DIVERGED_FNORM_NAN;
309189a9710SBarry Smith         PetscFunctionReturn(0);
310189a9710SBarry Smith       }
311534ebe21SPeter Brune     }
3123542a6bcSPeter Brune     /* Monitor convergence */
313ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
3143542a6bcSPeter Brune     snes->iter = i+1;
3153542a6bcSPeter Brune     snes->norm = fnorm;
316ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
317a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
3183542a6bcSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
3193542a6bcSPeter Brune     /* Test for convergence */
320365a6726SPeter Brune     if (normschedule == SNES_NORM_ALWAYS) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
3213542a6bcSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
3223542a6bcSPeter Brune     /* Call general purpose update function */
3233542a6bcSPeter Brune     if (snes->ops->update) {
3243542a6bcSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
3253542a6bcSPeter Brune     }
3263542a6bcSPeter Brune   }
327365a6726SPeter Brune   if (normschedule == SNES_NORM_ALWAYS) {
328534ebe21SPeter Brune     if (i == snes->max_its) {
329534ebe21SPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
330534ebe21SPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
331534ebe21SPeter Brune     }
3321aa26658SKarl Rupp   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
3333542a6bcSPeter Brune   PetscFunctionReturn(0);
3343542a6bcSPeter Brune }
3353542a6bcSPeter Brune 
3363542a6bcSPeter Brune /*MC
3373542a6bcSPeter Brune   SNESGS - Just calls the user-provided solution routine provided with SNESSetGS()
3383542a6bcSPeter Brune 
3393542a6bcSPeter Brune    Level: advanced
3403542a6bcSPeter Brune 
3413542a6bcSPeter Brune   Notes:
3423542a6bcSPeter Brune   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetPC(), it will have
3433542a6bcSPeter Brune   its parent's Gauss-Seidel routine associated with it.
3443542a6bcSPeter Brune 
3453542a6bcSPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types)
3463542a6bcSPeter Brune M*/
3473542a6bcSPeter Brune 
3483542a6bcSPeter Brune #undef __FUNCT__
3493542a6bcSPeter Brune #define __FUNCT__ "SNESCreate_GS"
3508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_GS(SNES snes)
3513542a6bcSPeter Brune {
3523542a6bcSPeter Brune   SNES_GS        *gs;
3533542a6bcSPeter Brune   PetscErrorCode ierr;
3543542a6bcSPeter Brune 
3553542a6bcSPeter Brune   PetscFunctionBegin;
3563542a6bcSPeter Brune   snes->ops->destroy        = SNESDestroy_GS;
3573542a6bcSPeter Brune   snes->ops->setup          = SNESSetUp_GS;
3583542a6bcSPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_GS;
3593542a6bcSPeter Brune   snes->ops->view           = SNESView_GS;
3603542a6bcSPeter Brune   snes->ops->solve          = SNESSolve_GS;
3613542a6bcSPeter Brune   snes->ops->reset          = SNESReset_GS;
3623542a6bcSPeter Brune 
3633542a6bcSPeter Brune   snes->usesksp = PETSC_FALSE;
3643542a6bcSPeter Brune   snes->usespc  = PETSC_FALSE;
3653542a6bcSPeter Brune 
36688976e71SPeter Brune   if (!snes->tolerancesset) {
3670e444f03SPeter Brune     snes->max_its   = 10000;
3680e444f03SPeter Brune     snes->max_funcs = 10000;
36988976e71SPeter Brune   }
3700e444f03SPeter Brune 
3713542a6bcSPeter Brune   ierr = PetscNewLog(snes, SNES_GS, &gs);CHKERRQ(ierr);
3726cd56b8bSPeter Brune 
3736cd56b8bSPeter Brune   gs->sweeps  = 1;
374b6266c6eSPeter Brune   gs->rtol    = 1e-5;
375b6266c6eSPeter Brune   gs->abstol  = 1e-15;
376b6266c6eSPeter Brune   gs->stol    = 1e-12;
377b6266c6eSPeter Brune   gs->max_its = 50;
378*737a7e12SPeter Brune   gs->h       = 1e-8;
3796cd56b8bSPeter Brune 
3803542a6bcSPeter Brune   snes->data = (void*) gs;
3813542a6bcSPeter Brune   PetscFunctionReturn(0);
3823542a6bcSPeter Brune }
383