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