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