xref: /petsc/src/snes/impls/gs/snesgs.c (revision f13dfd9ea68e0ddeee984e65c377a1819eab8a8a) !
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: [](ch_snes), `SNES`, `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: [](ch_snes), `SNES`, `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   Logically Collective
88 
89   Input Parameters:
90 + snes   - the `SNES` context
91 - sweeps - the number of sweeps of nonlinear GS to perform.
92 
93   Options Database Key:
94 . -snes_ngs_sweeps <n> - Number of sweeps of nonlinear GS to apply
95 
96   Level: intermediate
97 
98 .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSGetSweeps()`
99 @*/
100 PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
101 {
102   SNES_NGS *gs = (SNES_NGS *)snes->data;
103 
104   PetscFunctionBegin;
105   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
106   gs->sweeps = sweeps;
107   PetscFunctionReturn(PETSC_SUCCESS);
108 }
109 
110 /*@
111   SNESNGSGetSweeps - Gets the number of sweeps nonlinear GS will use in `SNESNCG`
112 
113   Input Parameter:
114 . snes - the `SNES` context
115 
116   Output Parameter:
117 . sweeps - the number of sweeps of nonlinear GS to perform.
118 
119   Level: intermediate
120 
121 .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESSetNGS()`, `SNESGetNGS()`, `SNESSetNPC()`, `SNESNGSSetSweeps()`
122 @*/
123 PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt *sweeps)
124 {
125   SNES_NGS *gs = (SNES_NGS *)snes->data;
126 
127   PetscFunctionBegin;
128   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
129   *sweeps = gs->sweeps;
130   PetscFunctionReturn(PETSC_SUCCESS);
131 }
132 
133 static PetscErrorCode SNESReset_NGS(SNES snes)
134 {
135   SNES_NGS *gs = (SNES_NGS *)snes->data;
136 
137   PetscFunctionBegin;
138   PetscCall(ISColoringDestroy(&gs->coloring));
139   PetscFunctionReturn(PETSC_SUCCESS);
140 }
141 
142 static PetscErrorCode SNESDestroy_NGS(SNES snes)
143 {
144   PetscFunctionBegin;
145   PetscCall(SNESReset_NGS(snes));
146   PetscCall(PetscFree(snes->data));
147   PetscFunctionReturn(PETSC_SUCCESS);
148 }
149 
150 static PetscErrorCode SNESSetUp_NGS(SNES snes)
151 {
152   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
153 
154   PetscFunctionBegin;
155   PetscCall(SNESGetNGS(snes, &f, NULL));
156   if (!f) PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
157   PetscFunctionReturn(PETSC_SUCCESS);
158 }
159 
160 static PetscErrorCode SNESSetFromOptions_NGS(SNES snes, PetscOptionItems *PetscOptionsObject)
161 {
162   SNES_NGS *gs = (SNES_NGS *)snes->data;
163   PetscInt  sweeps, max_its = PETSC_DEFAULT;
164   PetscReal rtol = PETSC_DEFAULT, atol = PETSC_DEFAULT, stol = PETSC_DEFAULT;
165   PetscBool flg, flg1, flg2, flg3;
166 
167   PetscFunctionBegin;
168   PetscOptionsHeadBegin(PetscOptionsObject, "SNES GS options");
169   /* GS Options */
170   PetscCall(PetscOptionsInt("-snes_ngs_sweeps", "Number of sweeps of GS to apply", "SNESComputeGS", gs->sweeps, &sweeps, &flg));
171   if (flg) PetscCall(SNESNGSSetSweeps(snes, sweeps));
172   PetscCall(PetscOptionsReal("-snes_ngs_atol", "Absolute residual tolerance for GS iteration", "SNESComputeGS", gs->abstol, &atol, &flg));
173   PetscCall(PetscOptionsReal("-snes_ngs_rtol", "Relative residual tolerance for GS iteration", "SNESComputeGS", gs->rtol, &rtol, &flg1));
174   PetscCall(PetscOptionsReal("-snes_ngs_stol", "Absolute update tolerance for GS iteration", "SNESComputeGS", gs->stol, &stol, &flg2));
175   PetscCall(PetscOptionsInt("-snes_ngs_max_it", "Maximum number of sweeps of GS to apply", "SNESComputeGS", gs->max_its, &max_its, &flg3));
176   if (flg || flg1 || flg2 || flg3) PetscCall(SNESNGSSetTolerances(snes, atol, rtol, stol, max_its));
177   flg = PETSC_FALSE;
178   PetscCall(PetscOptionsBool("-snes_ngs_secant", "Use finite difference secant approximation with coloring", "", flg, &flg, NULL));
179   if (flg) {
180     PetscCall(SNESSetNGS(snes, SNESComputeNGSDefaultSecant, NULL));
181     PetscCall(PetscInfo(snes, "Setting default finite difference secant approximation with coloring\n"));
182   }
183   PetscCall(PetscOptionsReal("-snes_ngs_secant_h", "Differencing parameter for secant search", "", gs->h, &gs->h, NULL));
184   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));
185 
186   PetscOptionsHeadEnd();
187   PetscFunctionReturn(PETSC_SUCCESS);
188 }
189 
190 static PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
191 {
192   PetscErrorCode (*f)(SNES, Vec, Vec, void *);
193   SNES_NGS *gs = (SNES_NGS *)snes->data;
194   PetscBool iascii;
195 
196   PetscFunctionBegin;
197   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
198   if (iascii) {
199     PetscCall(DMSNESGetNGS(snes->dm, &f, NULL));
200     if (f == SNESComputeNGSDefaultSecant) PetscCall(PetscViewerASCIIPrintf(viewer, "  Use finite difference secant approximation with coloring with h = %g \n", (double)gs->h));
201   }
202   PetscFunctionReturn(PETSC_SUCCESS);
203 }
204 
205 static PetscErrorCode SNESSolve_NGS(SNES snes)
206 {
207   Vec              F;
208   Vec              X;
209   Vec              B;
210   PetscInt         i;
211   PetscReal        fnorm;
212   SNESNormSchedule normschedule;
213 
214   PetscFunctionBegin;
215   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);
216 
217   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
218   X = snes->vec_sol;
219   F = snes->vec_func;
220   B = snes->vec_rhs;
221 
222   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
223   snes->iter = 0;
224   snes->norm = 0.;
225   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
226   snes->reason = SNES_CONVERGED_ITERATING;
227 
228   PetscCall(SNESGetNormSchedule(snes, &normschedule));
229   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
230     /* compute the initial function and preconditioned update delX */
231     if (!snes->vec_func_init_set) {
232       PetscCall(SNESComputeFunction(snes, X, F));
233     } else snes->vec_func_init_set = PETSC_FALSE;
234 
235     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
236     SNESCheckFunctionNorm(snes, fnorm);
237     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
238     snes->iter = 0;
239     snes->norm = fnorm;
240     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
241     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
242 
243     /* test convergence */
244     PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
245     PetscCall(SNESMonitor(snes, 0, snes->norm));
246     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
247   } else {
248     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
249     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
250   }
251 
252   /* Call general purpose update function */
253   PetscTryTypeMethod(snes, update, snes->iter);
254 
255   for (i = 0; i < snes->max_its; i++) {
256     PetscCall(SNESComputeNGS(snes, B, X));
257     /* only compute norms if requested or about to exit due to maximum iterations */
258     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
259       PetscCall(SNESComputeFunction(snes, X, F));
260       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
261       SNESCheckFunctionNorm(snes, fnorm);
262     }
263     /* Monitor convergence */
264     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
265     snes->iter = i + 1;
266     snes->norm = fnorm;
267     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
268     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
269     /* Test for convergence */
270     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
271     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
272     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
273     /* Call general purpose update function */
274     PetscTryTypeMethod(snes, update, snes->iter);
275   }
276   PetscFunctionReturn(PETSC_SUCCESS);
277 }
278 
279 /*MC
280   SNESNGS - Either calls the user-provided Gauss-Seidel solution routine provided with `SNESSetNGS()` or does a finite difference secant approximation
281             using coloring {cite}`bruneknepleysmithtu15`.
282 
283    Level: advanced
284 
285   Options Database Keys:
286 +   -snes_ngs_sweeps <n>                                                          - Number of sweeps of nonlinear GS to apply
287 .    -snes_ngs_atol <atol>                                                        - Absolute residual tolerance for nonlinear GS iteration
288 .    -snes_ngs_rtol <rtol>                                                        - Relative residual tolerance for nonlinear GS iteration
289 .    -snes_ngs_stol <stol>                                                        - Absolute update tolerance for nonlinear GS iteration
290 .    -snes_ngs_max_it <maxit>                                                     - Maximum number of sweeps of nonlinea GS to apply
291 .    -snes_ngs_secant                                                             - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine,
292                                                                                     this is used by default if no user provided Gauss-Seidel routine is available.
293                                                                                     Requires either that a `DM` that can compute a coloring
294                                                                                     is available or a Jacobian sparse matrix is provided (from which to get the coloring).
295 .    -snes_ngs_secant_h <h>                                                       - Differencing parameter for secant approximation
296 .    -snes_ngs_secant_mat_coloring                                                - Use the graph coloring of the Jacobian for the secant GS even if a `DM` is available.
297 -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed
298 
299   Notes:
300   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with `SNESGetNPC()`, it will have
301   its parent's Gauss-Seidel routine associated with it.
302 
303   By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with `SNESSetNormSchedule()`
304   or -snes_norm_schedule none
305 
306 .seealso: [](ch_snes), `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetNGS()`, `SNESType`, `SNESNGSSetSweeps()`, `SNESNGSSetTolerances()`,
307           `SNESSetNormSchedule()`, `SNESNGSGetTolerances()`, `SNESNGSSetSweeps()`
308 M*/
309 
310 PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
311 {
312   SNES_NGS *gs;
313 
314   PetscFunctionBegin;
315   snes->ops->destroy        = SNESDestroy_NGS;
316   snes->ops->setup          = SNESSetUp_NGS;
317   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
318   snes->ops->view           = SNESView_NGS;
319   snes->ops->solve          = SNESSolve_NGS;
320   snes->ops->reset          = SNESReset_NGS;
321 
322   snes->usesksp = PETSC_FALSE;
323   snes->usesnpc = PETSC_FALSE;
324 
325   snes->alwayscomputesfinalresidual = PETSC_FALSE;
326 
327   if (!snes->tolerancesset) {
328     snes->max_its   = 10000;
329     snes->max_funcs = 10000;
330   }
331 
332   PetscCall(PetscNew(&gs));
333 
334   gs->sweeps  = 1;
335   gs->rtol    = 1e-5;
336   gs->abstol  = PETSC_MACHINE_EPSILON;
337   gs->stol    = 1000 * PETSC_MACHINE_EPSILON;
338   gs->max_its = 50;
339   gs->h       = PETSC_SQRT_MACHINE_EPSILON;
340 
341   snes->data = (void *)gs;
342   PetscFunctionReturn(PETSC_SUCCESS);
343 }
344