xref: /petsc/src/snes/impls/gs/snesgs.c (revision 517f4e3460cd8c7e73c19b8bf533f1efe47f30a7)
1 #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h"  I*/
2 
3 /*@
4   SNESNGSSetTolerances - Sets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG`
5 
6   Logically Collective
7 
8   Input Parameters:
9 + snes   - the `SNES` context
10 . abstol - absolute convergence tolerance
11 . rtol   - relative convergence tolerance
12 . stol   - convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
13 - maxit  - maximum number of iterations
14 
15   Options Database Keys:
16 + -snes_ngs_atol <abstol> - Sets abstol
17 . -snes_ngs_rtol <rtol>   - Sets rtol
18 . -snes_ngs_stol <stol>   - Sets stol
19 - -snes_max_it <maxit>    - Sets maxit
20 
21   Level: intermediate
22 
23 .seealso: `SNESNCG`, `SNESSetTrustRegionTolerance()`
24 @*/
25 PetscErrorCode SNESNGSSetTolerances(SNES snes, PetscReal abstol, PetscReal rtol, PetscReal stol, PetscInt maxit)
26 {
27   SNES_NGS *gs = (SNES_NGS *)snes->data;
28 
29   PetscFunctionBegin;
30   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
31 
32   if (abstol != (PetscReal)PETSC_DEFAULT) {
33     PetscCheck(abstol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %g must be non-negative", (double)abstol);
34     gs->abstol = abstol;
35   }
36   if (rtol != (PetscReal)PETSC_DEFAULT) {
37     PetscCheck(rtol >= 0.0 && rtol < 1.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %g must be non-negative and less than 1.0", (double)rtol);
38     gs->rtol = rtol;
39   }
40   if (stol != (PetscReal)PETSC_DEFAULT) {
41     PetscCheck(stol >= 0.0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Step tolerance %g must be non-negative", (double)stol);
42     gs->stol = stol;
43   }
44   if (maxit != PETSC_DEFAULT) {
45     PetscCheck(maxit >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", maxit);
46     gs->max_its = maxit;
47   }
48   PetscFunctionReturn(PETSC_SUCCESS);
49 }
50 
51 /*@
52   SNESNGSGetTolerances - Gets various parameters used in convergence tests for nonlinear Gauss-Seidel `SNESNCG`
53 
54   Not Collective
55 
56   Input Parameters:
57 + snes  - the `SNES` context
58 . atol  - absolute convergence tolerance
59 . rtol  - relative convergence tolerance
60 . stol  - convergence tolerance in terms of the norm
61            of the change in the solution between steps
62 - maxit - maximum number of iterations
63 
64   Level: intermediate
65 
66   Note:
67   The user can specify NULL for any parameter that is not needed.
68 
69 .seealso: `SNESNCG`, `SNESSetTolerances()`
70 @*/
71 PetscErrorCode SNESNGSGetTolerances(SNES snes, PetscReal *atol, PetscReal *rtol, PetscReal *stol, PetscInt *maxit)
72 {
73   SNES_NGS *gs = (SNES_NGS *)snes->data;
74 
75   PetscFunctionBegin;
76   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
77   if (atol) *atol = gs->abstol;
78   if (rtol) *rtol = gs->rtol;
79   if (stol) *stol = gs->stol;
80   if (maxit) *maxit = gs->max_its;
81   PetscFunctionReturn(PETSC_SUCCESS);
82 }
83 
84 /*@
85   SNESNGSSetSweeps - Sets the number of sweeps of nonlinear GS to use in `SNESNCG`
86 
87   Input Parameters:
88 + snes   - the `SNES` context
89 - sweeps - the number of sweeps of nonlinear GS to perform.
90 
91   Options Database Key:
92 . -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply
93 
94   Level: intermediate
95 
96 .seealso: `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()`
97 @*/
98 PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
99 {
100   SNES_NGS *gs = (SNES_NGS *)snes->data;
101 
102   PetscFunctionBegin;
103   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
104   gs->sweeps = sweeps;
105   PetscFunctionReturn(PETSC_SUCCESS);
106 }
107 
108 /*@
109   SNESNGSGetSweeps - Gets the number of sweeps nonlinear GS will use in `SNESNCG`
110 
111   Input Parameter:
112 . snes - the `SNES` context
113 
114   Output Parameter:
115 . sweeps - the number of sweeps of nonlinear GS to perform.
116 
117   Level: intermediate
118 
119 .seealso: `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()`
120 @*/
121 PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps)
122 {
123   SNES_NGS *gs = (SNES_NGS *)snes->data;
124 
125   PetscFunctionBegin;
126   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
127   *sweeps = gs->sweeps;
128   PetscFunctionReturn(PETSC_SUCCESS);
129 }
130 
131 static PetscErrorCode SNESReset_NGS(SNES snes)
132 {
133   SNES_NGS *gs = (SNES_NGS *)snes->data;
134 
135   PetscFunctionBegin;
136   PetscCall(ISColoringDestroy(&gs->coloring));
137   PetscFunctionReturn(PETSC_SUCCESS);
138 }
139 
140 static PetscErrorCode SNESDestroy_NGS(SNES snes)
141 {
142   PetscFunctionBegin;
143   PetscCall(SNESReset_NGS(snes));
144   PetscCall(PetscFree(snes->data));
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147 
148 static PetscErrorCode SNESSetUp_NGS(SNES snes)
149 {
150   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
151 
152   PetscFunctionBegin;
153   PetscCall(SNESGetNGS(snes, &f, NULL));
154   if (!f) PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
155   PetscFunctionReturn(PETSC_SUCCESS);
156 }
157 
158 static PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems *PetscOptionsObject)
159 {
160   SNES_NGS *gs = (SNES_NGS *)snes->data;
161   PetscInt  sweeps, max_its = PETSC_DEFAULT;
162   PetscReal rtol = PETSC_DEFAULT, atol = PETSC_DEFAULT, stol = PETSC_DEFAULT;
163   PetscBool flg, flg1, flg2, flg3;
164 
165   PetscFunctionBegin;
166   PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options");
167   /* GS Options */
168   PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg));
169   if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps));
170   PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg));
171   PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1));
172   PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2));
173   PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3));
174   if (flg || flg1 || flg2 || flg3) PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its));
175   flg = PETSC_FALSE;
176   PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL));
177   if (flg) {
178     PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
179     PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n"));
180   }
181   PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL));
182   PetscCall(PetscOptionsBool("-snes_ngs_secant_mat_coloring", "Use the graph coloring of the Jacobian for the secant GS", "", gs->secant_mat, &gs->secant_mat, &flg));
183 
184   PetscOptionsHeadEnd();
185   PetscFunctionReturn(PETSC_SUCCESS);
186 }
187 
188 static PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
189 {
190   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
191   SNES_NGS *gs = (SNES_NGS *)snes->data;
192   PetscBool iascii;
193 
194   PetscFunctionBegin;
195   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
196   if (iascii) {
197     PetscCall(DMSNESGetNGS(snes->dm, &f, NULL));
198     if (f == SNESComputeNGSDefaultSecant) PetscCall(PetscViewerASCIIPrintf(viewer, "  Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h));
199   }
200   PetscFunctionReturn(PETSC_SUCCESS);
201 }
202 
203 static PetscErrorCode SNESSolve_NGS(SNES snes)
204 {
205   Vec              F;
206   Vec              X;
207   Vec              B;
208   PetscInt         i;
209   PetscReal        fnorm;
210   SNESNormSchedule normschedule;
211 
212   PetscFunctionBegin;
213 
214   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
215 
216   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
217   X = snes->vec_sol;
218   F = snes->vec_func;
219   B = snes->vec_rhs;
220 
221   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
222   snes->iter = 0;
223   snes->norm = 0.;
224   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
225   snes->reason = SNES_CONVERGED_ITERATING;
226 
227   PetscCall(SNESGetNormSchedule(snes, &normschedule));
228   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
229     /* compute the initial function and preconditioned update delX */
230     if (!snes->vec_func_init_set) {
231       PetscCall(SNESComputeFunction(snes, X, F));
232     } else snes->vec_func_init_set = PETSC_FALSE;
233 
234     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
235     SNESCheckFunctionNorm(snes, fnorm);
236     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
237     snes->iter = 0;
238     snes->norm = fnorm;
239     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
240     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
241 
242     /* test convergence */
243     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
244     PetscCall(SNESMonitor(snes, 0, snes->norm));
245     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
246   } else {
247     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
248     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
249   }
250 
251   /* Call general purpose update function */
252   PetscTryTypeMethod(snes, update, snes->iter);
253 
254   for (i = 0; i < snes->max_its; i++) {
255     PetscCall(SNESComputeNGS(snes, B, X));
256     /* only compute norms if requested or about to exit due to maximum iterations */
257     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
258       PetscCall(SNESComputeFunction(snes, X, F));
259       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
260       SNESCheckFunctionNorm(snes, fnorm);
261     }
262     /* Monitor convergence */
263     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
264     snes->iter = i + 1;
265     snes->norm = fnorm;
266     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
267     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
268     /* Test for convergence */
269     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
270     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
271     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
272     /* Call general purpose update function */
273     PetscTryTypeMethod(snes, update, snes->iter);
274   }
275   PetscFunctionReturn(PETSC_SUCCESS);
276 }
277 
278 /*MC
279   SNESNGS - Either calls the user-provided solution routine provided with `SNESSetNGS()` or does a finite difference secant approximation
280             using coloring.
281 
282    Level: advanced
283 
284   Options Database Keys:
285 +   -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply
286 .    -snes_ngs_atol <atol> - Absolute residual tolerance for nonlinear GS iteration
287 .    -snes_ngs_rtol <rtol> - Relative residual tolerance for nonlinear GS iteration
288 .    -snes_ngs_stol <stol> - Absolute update tolerance for nonlinear GS iteration
289 .    -snes_ngs_max_it <maxit> - Maximum number of sweeps of nonlinea GS to apply
290 .    -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is
291                         used by default if no user provided Gauss-Seidel routine is available. Requires either that a `DM` that can compute a coloring
292                         is available or a Jacobian sparse matrix is provided (from which to get the coloring).
293 .    -snes_ngs_secant_h <h> - Differencing parameter for secant approximation
294 .    -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available.
295 -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed
296 
297   Notes:
298   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with `SNESGetNPC()`, it will have
299   its parent's Gauss-Seidel routine associated with it.
300 
301   By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with `SNESSetNormSchedule()`
302   or -snes_norm_schedule none
303 
304    References:
305 .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
306    SIAM Review, 57(4), 2015
307 
308 .seealso: `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`,
309           `SNESSetNormSchedule()`
310 M*/
311 
312 PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
313 {
314   SNES_NGS *gs;
315 
316   PetscFunctionBegin;
317   snes->ops->destroy        = SNESDestroy_NGS;
318   snes->ops->setup          = SNESSetUp_NGS;
319   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
320   snes->ops->view           = SNESView_NGS;
321   snes->ops->solve          = SNESSolve_NGS;
322   snes->ops->reset          = SNESReset_NGS;
323 
324   snes->usesksp = PETSC_FALSE;
325   snes->usesnpc = PETSC_FALSE;
326 
327   snes->alwayscomputesfinalresidual = PETSC_FALSE;
328 
329   if (!snes->tolerancesset) {
330     snes->max_its   = 10000;
331     snes->max_funcs = 10000;
332   }
333 
334   PetscCall(PetscNew(&gs));
335 
336   gs->sweeps  = 1;
337   gs->rtol    = 1e-5;
338   gs->abstol  = PETSC_MACHINE_EPSILON;
339   gs->stol    = 1000 * PETSC_MACHINE_EPSILON;
340   gs->max_its = 50;
341   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
342 
343   snes->data = (void *)gs;
344   PetscFunctionReturn(PETSC_SUCCESS);
345 }
346