xref: /petsc/src/snes/interface/snes.c (revision c5589a102014197226e080e0d4371932bf22410e)
1 /*$Id: snes.c,v 1.235 2001/08/21 21:03:49 bsmith Exp $*/
2 
3 #include "src/snes/snesimpl.h"      /*I "petscsnes.h"  I*/
4 
5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6 PetscFList SNESList              = PETSC_NULL;
7 
8 /* Logging support */
9 int SNES_COOKIE = 0;
10 int SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0;
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "SNESView"
14 /*@C
15    SNESView - Prints the SNES data structure.
16 
17    Collective on SNES
18 
19    Input Parameters:
20 +  SNES - the SNES context
21 -  viewer - visualization context
22 
23    Options Database Key:
24 .  -snes_view - Calls SNESView() at end of SNESSolve()
25 
26    Notes:
27    The available visualization contexts include
28 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
29 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
30          output where only the first processor opens
31          the file.  All other processors send their
32          data to the first processor to print.
33 
34    The user can open an alternative visualization context with
35    PetscViewerASCIIOpen() - output to a specified file.
36 
37    Level: beginner
38 
39 .keywords: SNES, view
40 
41 .seealso: PetscViewerASCIIOpen()
42 @*/
43 int SNESView(SNES snes,PetscViewer viewer)
44 {
45   SNES_KSP_EW_ConvCtx *kctx;
46   int                 ierr;
47   SLES                sles;
48   char                *type;
49   PetscTruth          isascii,isstring;
50 
51   PetscFunctionBegin;
52   PetscValidHeaderSpecific(snes,SNES_COOKIE);
53   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);
54   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
55   PetscCheckSameComm(snes,viewer);
56 
57   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
58   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
59   if (isascii) {
60     if (snes->prefix) {
61       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr);
62     } else {
63       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
64     }
65     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
66     if (type) {
67       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
68     } else {
69       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
70     }
71     if (snes->view) {
72       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
73       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
74       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
75     }
76     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
77     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",
78                  snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr);
79     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
80     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
81     if (snes->ksp_ewconv) {
82       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
83       if (kctx) {
84         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
85         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
86         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
87       }
88     }
89   } else if (isstring) {
90     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
91     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
92   }
93   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
94   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
95   ierr = SLESView(sles,viewer);CHKERRQ(ierr);
96   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 /*
101   We retain a list of functions that also take SNES command
102   line options. These are called at the end SNESSetFromOptions()
103 */
104 #define MAXSETFROMOPTIONS 5
105 static int numberofsetfromoptions;
106 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "SNESAddOptionsChecker"
110 /*@
111   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
112 
113   Not Collective
114 
115   Input Parameter:
116 . snescheck - function that checks for options
117 
118   Level: developer
119 
120 .seealso: SNESSetFromOptions()
121 @*/
122 int SNESAddOptionsChecker(int (*snescheck)(SNES))
123 {
124   PetscFunctionBegin;
125   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
126     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS);
127   }
128 
129   othersetfromoptions[numberofsetfromoptions++] = snescheck;
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "SNESSetFromOptions"
135 /*@
136    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
137 
138    Collective on SNES
139 
140    Input Parameter:
141 .  snes - the SNES context
142 
143    Options Database Keys:
144 +  -snes_type <type> - ls, tr, umls, umtr, test
145 .  -snes_stol - convergence tolerance in terms of the norm
146                 of the change in the solution between steps
147 .  -snes_atol <atol> - absolute tolerance of residual norm
148 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
149 .  -snes_max_it <max_it> - maximum number of iterations
150 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
151 .  -snes_max_fail <max_fail> - maximum number of failures
152 .  -snes_trtol <trtol> - trust region tolerance
153 .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
154                                solver; hence iterations will continue until max_it
155                                or some other criterion is reached. Saves expense
156                                of convergence test
157 .  -snes_monitor - prints residual norm at each iteration
158 .  -snes_vecmonitor - plots solution at each iteration
159 .  -snes_vecmonitor_update - plots update to solution at each iteration
160 .  -snes_xmonitor - plots residual norm at each iteration
161 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
162 -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
163 
164     Options Database for Eisenstat-Walker method:
165 +  -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence
166 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
167 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
168 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
169 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
170 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
171 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
172 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
173 
174    Notes:
175    To see all options, run your program with the -help option or consult
176    the users manual.
177 
178    Level: beginner
179 
180 .keywords: SNES, nonlinear, set, options, database
181 
182 .seealso: SNESSetOptionsPrefix()
183 @*/
184 int SNESSetFromOptions(SNES snes)
185 {
186   SLES                sles;
187   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
188   PetscTruth          flg;
189   int                 ierr, i;
190   char                *deft,type[256];
191 
192   PetscFunctionBegin;
193   PetscValidHeaderSpecific(snes,SNES_COOKIE);
194 
195   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
196     if (snes->type_name) {
197       deft = snes->type_name;
198     } else {
199       deft = SNESLS;
200     }
201 
202     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
203     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
204     if (flg) {
205       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
206     } else if (!snes->type_name) {
207       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
208     }
209     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
210 
211     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
212     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
213 
214     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
215     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
216     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
217     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
218 
219     ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
220 
221     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
222     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
223     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
224     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
225     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
226     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
227     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
228 
229     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
230     if (flg) {snes->converged = 0;}
231     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
232     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
233     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
234     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
235     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
236     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
237     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
238     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
239     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
240     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
241     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
242     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
243     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
244     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
245     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
246     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
247 
248     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
249     if (flg) {
250       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
251       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
252     }
253 
254     for(i = 0; i < numberofsetfromoptions; i++) {
255       ierr = (*othersetfromoptions[i])(snes);                                                             CHKERRQ(ierr);
256     }
257 
258     if (snes->setfromoptions) {
259       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
260     }
261 
262   ierr = PetscOptionsEnd();CHKERRQ(ierr);
263 
264   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
265   ierr = SLESSetFromOptions(sles);CHKERRQ(ierr);
266 
267   PetscFunctionReturn(0);
268 }
269 
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "SNESSetApplicationContext"
273 /*@
274    SNESSetApplicationContext - Sets the optional user-defined context for
275    the nonlinear solvers.
276 
277    Collective on SNES
278 
279    Input Parameters:
280 +  snes - the SNES context
281 -  usrP - optional user context
282 
283    Level: intermediate
284 
285 .keywords: SNES, nonlinear, set, application, context
286 
287 .seealso: SNESGetApplicationContext()
288 @*/
289 int SNESSetApplicationContext(SNES snes,void *usrP)
290 {
291   PetscFunctionBegin;
292   PetscValidHeaderSpecific(snes,SNES_COOKIE);
293   snes->user		= usrP;
294   PetscFunctionReturn(0);
295 }
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "SNESGetApplicationContext"
299 /*@C
300    SNESGetApplicationContext - Gets the user-defined context for the
301    nonlinear solvers.
302 
303    Not Collective
304 
305    Input Parameter:
306 .  snes - SNES context
307 
308    Output Parameter:
309 .  usrP - user context
310 
311    Level: intermediate
312 
313 .keywords: SNES, nonlinear, get, application, context
314 
315 .seealso: SNESSetApplicationContext()
316 @*/
317 int SNESGetApplicationContext(SNES snes,void **usrP)
318 {
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(snes,SNES_COOKIE);
321   *usrP = snes->user;
322   PetscFunctionReturn(0);
323 }
324 
325 #undef __FUNCT__
326 #define __FUNCT__ "SNESGetIterationNumber"
327 /*@
328    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
329    at this time.
330 
331    Not Collective
332 
333    Input Parameter:
334 .  snes - SNES context
335 
336    Output Parameter:
337 .  iter - iteration number
338 
339    Notes:
340    For example, during the computation of iteration 2 this would return 1.
341 
342    This is useful for using lagged Jacobians (where one does not recompute the
343    Jacobian at each SNES iteration). For example, the code
344 .vb
345       ierr = SNESGetIterationNumber(snes,&it);
346       if (!(it % 2)) {
347         [compute Jacobian here]
348       }
349 .ve
350    can be used in your ComputeJacobian() function to cause the Jacobian to be
351    recomputed every second SNES iteration.
352 
353    Level: intermediate
354 
355 .keywords: SNES, nonlinear, get, iteration, number
356 @*/
357 int SNESGetIterationNumber(SNES snes,int* iter)
358 {
359   PetscFunctionBegin;
360   PetscValidHeaderSpecific(snes,SNES_COOKIE);
361   PetscValidIntPointer(iter);
362   *iter = snes->iter;
363   PetscFunctionReturn(0);
364 }
365 
366 #undef __FUNCT__
367 #define __FUNCT__ "SNESGetFunctionNorm"
368 /*@
369    SNESGetFunctionNorm - Gets the norm of the current function that was set
370    with SNESSSetFunction().
371 
372    Collective on SNES
373 
374    Input Parameter:
375 .  snes - SNES context
376 
377    Output Parameter:
378 .  fnorm - 2-norm of function
379 
380    Level: intermediate
381 
382 .keywords: SNES, nonlinear, get, function, norm
383 
384 .seealso: SNESGetFunction()
385 @*/
386 int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
387 {
388   PetscFunctionBegin;
389   PetscValidHeaderSpecific(snes,SNES_COOKIE);
390   PetscValidScalarPointer(fnorm);
391   *fnorm = snes->norm;
392   PetscFunctionReturn(0);
393 }
394 
395 #undef __FUNCT__
396 #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
397 /*@
398    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
399    attempted by the nonlinear solver.
400 
401    Not Collective
402 
403    Input Parameter:
404 .  snes - SNES context
405 
406    Output Parameter:
407 .  nfails - number of unsuccessful steps attempted
408 
409    Notes:
410    This counter is reset to zero for each successive call to SNESSolve().
411 
412    Level: intermediate
413 
414 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
415 @*/
416 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
417 {
418   PetscFunctionBegin;
419   PetscValidHeaderSpecific(snes,SNES_COOKIE);
420   PetscValidIntPointer(nfails);
421   *nfails = snes->numFailures;
422   PetscFunctionReturn(0);
423 }
424 
425 #undef __FUNCT__
426 #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
427 /*@
428    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
429    attempted by the nonlinear solver before it gives up.
430 
431    Not Collective
432 
433    Input Parameters:
434 +  snes     - SNES context
435 -  maxFails - maximum of unsuccessful steps
436 
437    Level: intermediate
438 
439 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
440 @*/
441 int SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails)
442 {
443   PetscFunctionBegin;
444   PetscValidHeaderSpecific(snes,SNES_COOKIE);
445   snes->maxFailures = maxFails;
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
451 /*@
452    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
453    attempted by the nonlinear solver before it gives up.
454 
455    Not Collective
456 
457    Input Parameter:
458 .  snes     - SNES context
459 
460    Output Parameter:
461 .  maxFails - maximum of unsuccessful steps
462 
463    Level: intermediate
464 
465 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
466 @*/
467 int SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails)
468 {
469   PetscFunctionBegin;
470   PetscValidHeaderSpecific(snes,SNES_COOKIE);
471   PetscValidIntPointer(maxFails);
472   *maxFails = snes->maxFailures;
473   PetscFunctionReturn(0);
474 }
475 
476 #undef __FUNCT__
477 #define __FUNCT__ "SNESGetNumberLinearIterations"
478 /*@
479    SNESGetNumberLinearIterations - Gets the total number of linear iterations
480    used by the nonlinear solver.
481 
482    Not Collective
483 
484    Input Parameter:
485 .  snes - SNES context
486 
487    Output Parameter:
488 .  lits - number of linear iterations
489 
490    Notes:
491    This counter is reset to zero for each successive call to SNESSolve().
492 
493    Level: intermediate
494 
495 .keywords: SNES, nonlinear, get, number, linear, iterations
496 @*/
497 int SNESGetNumberLinearIterations(SNES snes,int* lits)
498 {
499   PetscFunctionBegin;
500   PetscValidHeaderSpecific(snes,SNES_COOKIE);
501   PetscValidIntPointer(lits);
502   *lits = snes->linear_its;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "SNESGetSLES"
508 /*@C
509    SNESGetSLES - Returns the SLES context for a SNES solver.
510 
511    Not Collective, but if SNES object is parallel, then SLES object is parallel
512 
513    Input Parameter:
514 .  snes - the SNES context
515 
516    Output Parameter:
517 .  sles - the SLES context
518 
519    Notes:
520    The user can then directly manipulate the SLES context to set various
521    options, etc.  Likewise, the user can then extract and manipulate the
522    KSP and PC contexts as well.
523 
524    Level: beginner
525 
526 .keywords: SNES, nonlinear, get, SLES, context
527 
528 .seealso: KSPGetPC(), SLESGetKSP()
529 @*/
530 int SNESGetSLES(SNES snes,SLES *sles)
531 {
532   PetscFunctionBegin;
533   PetscValidHeaderSpecific(snes,SNES_COOKIE);
534   *sles = snes->sles;
535   PetscFunctionReturn(0);
536 }
537 
538 #undef __FUNCT__
539 #define __FUNCT__ "SNESPublish_Petsc"
540 static int SNESPublish_Petsc(PetscObject obj)
541 {
542 #if defined(PETSC_HAVE_AMS)
543   SNES          v = (SNES) obj;
544   int          ierr;
545 #endif
546 
547   PetscFunctionBegin;
548 
549 #if defined(PETSC_HAVE_AMS)
550   /* if it is already published then return */
551   if (v->amem >=0) PetscFunctionReturn(0);
552 
553   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
554   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
555                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
556   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
557                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
558   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
559 #endif
560   PetscFunctionReturn(0);
561 }
562 
563 /* -----------------------------------------------------------*/
564 #undef __FUNCT__
565 #define __FUNCT__ "SNESCreate"
566 /*@C
567    SNESCreate - Creates a nonlinear solver context.
568 
569    Collective on MPI_Comm
570 
571    Input Parameters:
572 +  comm - MPI communicator
573 
574    Output Parameter:
575 .  outsnes - the new SNES context
576 
577    Options Database Keys:
578 +   -snes_mf - Activates default matrix-free Jacobian-vector products,
579                and no preconditioning matrix
580 .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
581                products, and a user-provided preconditioning matrix
582                as set by SNESSetJacobian()
583 -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
584 
585    Level: beginner
586 
587 .keywords: SNES, nonlinear, create, context
588 
589 .seealso: SNESSolve(), SNESDestroy(), SNES
590 @*/
591 int SNESCreate(MPI_Comm comm,SNES *outsnes)
592 {
593   int                 ierr;
594   SNES                snes;
595   SNES_KSP_EW_ConvCtx *kctx;
596 
597   PetscFunctionBegin;
598   PetscValidPointer(outsnes);
599   *outsnes = PETSC_NULL;
600 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
601   ierr = SNESInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
602 #endif
603 
604   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
605   PetscLogObjectCreate(snes);
606   snes->bops->publish     = SNESPublish_Petsc;
607   snes->max_its           = 50;
608   snes->max_funcs	  = 10000;
609   snes->norm		  = 0.0;
610   snes->rtol		  = 1.e-8;
611   snes->ttol              = 0.0;
612   snes->atol		  = 1.e-50;
613   snes->xtol		  = 1.e-8;
614   snes->deltatol	  = 1.e-12;
615   snes->nfuncs            = 0;
616   snes->numFailures       = 0;
617   snes->maxFailures       = 1;
618   snes->linear_its        = 0;
619   snes->numbermonitors    = 0;
620   snes->data              = 0;
621   snes->view              = 0;
622   snes->setupcalled       = 0;
623   snes->ksp_ewconv        = PETSC_FALSE;
624   snes->vwork             = 0;
625   snes->nwork             = 0;
626   snes->conv_hist_len     = 0;
627   snes->conv_hist_max     = 0;
628   snes->conv_hist         = PETSC_NULL;
629   snes->conv_hist_its     = PETSC_NULL;
630   snes->conv_hist_reset   = PETSC_TRUE;
631   snes->reason            = SNES_CONVERGED_ITERATING;
632 
633   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
634   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
635   PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
636   snes->kspconvctx  = (void*)kctx;
637   kctx->version     = 2;
638   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
639                              this was too large for some test cases */
640   kctx->rtol_last   = 0;
641   kctx->rtol_max    = .9;
642   kctx->gamma       = 1.0;
643   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
644   kctx->alpha       = kctx->alpha2;
645   kctx->threshold   = .1;
646   kctx->lresid_last = 0;
647   kctx->norm_last   = 0;
648 
649   ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr);
650   PetscLogObjectParent(snes,snes->sles)
651 
652   *outsnes = snes;
653   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
654   PetscFunctionReturn(0);
655 }
656 
657 /* --------------------------------------------------------------- */
658 #undef __FUNCT__
659 #define __FUNCT__ "SNESSetFunction"
660 /*@C
661    SNESSetFunction - Sets the function evaluation routine and function
662    vector for use by the SNES routines in solving systems of nonlinear
663    equations.
664 
665    Collective on SNES
666 
667    Input Parameters:
668 +  snes - the SNES context
669 .  func - function evaluation routine
670 .  r - vector to store function value
671 -  ctx - [optional] user-defined context for private data for the
672          function evaluation routine (may be PETSC_NULL)
673 
674    Calling sequence of func:
675 $    func (SNES snes,Vec x,Vec f,void *ctx);
676 
677 .  f - function vector
678 -  ctx - optional user-defined function context
679 
680    Notes:
681    The Newton-like methods typically solve linear systems of the form
682 $      f'(x) x = -f(x),
683    where f'(x) denotes the Jacobian matrix and f(x) is the function.
684 
685    Level: beginner
686 
687 .keywords: SNES, nonlinear, set, function
688 
689 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
690 @*/
691 int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
692 {
693   PetscFunctionBegin;
694   PetscValidHeaderSpecific(snes,SNES_COOKIE);
695   PetscValidHeaderSpecific(r,VEC_COOKIE);
696   PetscCheckSameComm(snes,r);
697 
698   snes->computefunction     = func;
699   snes->vec_func            = snes->vec_func_always = r;
700   snes->funP                = ctx;
701   PetscFunctionReturn(0);
702 }
703 
704 #undef __FUNCT__
705 #define __FUNCT__ "SNESComputeFunction"
706 /*@
707    SNESComputeFunction - Calls the function that has been set with
708                          SNESSetFunction().
709 
710    Collective on SNES
711 
712    Input Parameters:
713 +  snes - the SNES context
714 -  x - input vector
715 
716    Output Parameter:
717 .  y - function vector, as set by SNESSetFunction()
718 
719    Notes:
720    SNESComputeFunction() is typically used within nonlinear solvers
721    implementations, so most users would not generally call this routine
722    themselves.
723 
724    Level: developer
725 
726 .keywords: SNES, nonlinear, compute, function
727 
728 .seealso: SNESSetFunction(), SNESGetFunction()
729 @*/
730 int SNESComputeFunction(SNES snes,Vec x,Vec y)
731 {
732   int    ierr;
733 
734   PetscFunctionBegin;
735   PetscValidHeaderSpecific(snes,SNES_COOKIE);
736   PetscValidHeaderSpecific(x,VEC_COOKIE);
737   PetscValidHeaderSpecific(y,VEC_COOKIE);
738   PetscCheckSameComm(snes,x);
739   PetscCheckSameComm(snes,y);
740 
741   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
742   PetscStackPush("SNES user function");
743   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
744   PetscStackPop;
745   snes->nfuncs++;
746   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
747   PetscFunctionReturn(0);
748 }
749 
750 #undef __FUNCT__
751 #define __FUNCT__ "SNESComputeJacobian"
752 /*@
753    SNESComputeJacobian - Computes the Jacobian matrix that has been
754    set with SNESSetJacobian().
755 
756    Collective on SNES and Mat
757 
758    Input Parameters:
759 +  snes - the SNES context
760 -  x - input vector
761 
762    Output Parameters:
763 +  A - Jacobian matrix
764 .  B - optional preconditioning matrix
765 -  flag - flag indicating matrix structure
766 
767    Notes:
768    Most users should not need to explicitly call this routine, as it
769    is used internally within the nonlinear solvers.
770 
771    See SLESSetOperators() for important information about setting the
772    flag parameter.
773 
774    Level: developer
775 
776 .keywords: SNES, compute, Jacobian, matrix
777 
778 .seealso:  SNESSetJacobian(), SLESSetOperators()
779 @*/
780 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
781 {
782   int    ierr;
783 
784   PetscFunctionBegin;
785   PetscValidHeaderSpecific(snes,SNES_COOKIE);
786   PetscValidHeaderSpecific(X,VEC_COOKIE);
787   PetscCheckSameComm(snes,X);
788   if (!snes->computejacobian) PetscFunctionReturn(0);
789   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
790   *flg = DIFFERENT_NONZERO_PATTERN;
791   PetscStackPush("SNES user Jacobian function");
792   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
793   PetscStackPop;
794   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
795   /* make sure user returned a correct Jacobian and preconditioner */
796   PetscValidHeaderSpecific(*A,MAT_COOKIE);
797   PetscValidHeaderSpecific(*B,MAT_COOKIE);
798   PetscFunctionReturn(0);
799 }
800 
801 #undef __FUNCT__
802 #define __FUNCT__ "SNESSetJacobian"
803 /*@C
804    SNESSetJacobian - Sets the function to compute Jacobian as well as the
805    location to store the matrix.
806 
807    Collective on SNES and Mat
808 
809    Input Parameters:
810 +  snes - the SNES context
811 .  A - Jacobian matrix
812 .  B - preconditioner matrix (usually same as the Jacobian)
813 .  func - Jacobian evaluation routine
814 -  ctx - [optional] user-defined context for private data for the
815          Jacobian evaluation routine (may be PETSC_NULL)
816 
817    Calling sequence of func:
818 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
819 
820 +  x - input vector
821 .  A - Jacobian matrix
822 .  B - preconditioner matrix, usually the same as A
823 .  flag - flag indicating information about the preconditioner matrix
824    structure (same as flag in SLESSetOperators())
825 -  ctx - [optional] user-defined Jacobian context
826 
827    Notes:
828    See SLESSetOperators() for important information about setting the flag
829    output parameter in the routine func().  Be sure to read this information!
830 
831    The routine func() takes Mat * as the matrix arguments rather than Mat.
832    This allows the Jacobian evaluation routine to replace A and/or B with a
833    completely new new matrix structure (not just different matrix elements)
834    when appropriate, for instance, if the nonzero structure is changing
835    throughout the global iterations.
836 
837    Level: beginner
838 
839 .keywords: SNES, nonlinear, set, Jacobian, matrix
840 
841 .seealso: SLESSetOperators(), SNESSetFunction()
842 @*/
843 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
844 {
845   int ierr;
846 
847   PetscFunctionBegin;
848   PetscValidHeaderSpecific(snes,SNES_COOKIE);
849   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE);
850   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE);
851   if (A) PetscCheckSameComm(snes,A);
852   if (B) PetscCheckSameComm(snes,B);
853   if (func) snes->computejacobian = func;
854   if (ctx)  snes->jacP            = ctx;
855   if (A) {
856     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
857     snes->jacobian = A;
858     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
859   }
860   if (B) {
861     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
862     snes->jacobian_pre = B;
863     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
864   }
865   PetscFunctionReturn(0);
866 }
867 
868 #undef __FUNCT__
869 #define __FUNCT__ "SNESGetJacobian"
870 /*@C
871    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
872    provided context for evaluating the Jacobian.
873 
874    Not Collective, but Mat object will be parallel if SNES object is
875 
876    Input Parameter:
877 .  snes - the nonlinear solver context
878 
879    Output Parameters:
880 +  A - location to stash Jacobian matrix (or PETSC_NULL)
881 .  B - location to stash preconditioner matrix (or PETSC_NULL)
882 .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
883 -  func - location to put Jacobian function (or PETSC_NULL)
884 
885    Level: advanced
886 
887 .seealso: SNESSetJacobian(), SNESComputeJacobian()
888 @*/
889 int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
890 {
891   PetscFunctionBegin;
892   PetscValidHeaderSpecific(snes,SNES_COOKIE);
893   if (A)    *A    = snes->jacobian;
894   if (B)    *B    = snes->jacobian_pre;
895   if (ctx)  *ctx  = snes->jacP;
896   if (func) *func = snes->computejacobian;
897   PetscFunctionReturn(0);
898 }
899 
900 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
901 extern int SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
902 
903 #undef __FUNCT__
904 #define __FUNCT__ "SNESSetUp"
905 /*@
906    SNESSetUp - Sets up the internal data structures for the later use
907    of a nonlinear solver.
908 
909    Collective on SNES
910 
911    Input Parameters:
912 +  snes - the SNES context
913 -  x - the solution vector
914 
915    Notes:
916    For basic use of the SNES solvers the user need not explicitly call
917    SNESSetUp(), since these actions will automatically occur during
918    the call to SNESSolve().  However, if one wishes to control this
919    phase separately, SNESSetUp() should be called after SNESCreate()
920    and optional routines of the form SNESSetXXX(), but before SNESSolve().
921 
922    Level: advanced
923 
924 .keywords: SNES, nonlinear, setup
925 
926 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
927 @*/
928 int SNESSetUp(SNES snes,Vec x)
929 {
930   int        ierr;
931   PetscTruth flg, iseqtr;
932 
933   PetscFunctionBegin;
934   PetscValidHeaderSpecific(snes,SNES_COOKIE);
935   PetscValidHeaderSpecific(x,VEC_COOKIE);
936   PetscCheckSameComm(snes,x);
937   snes->vec_sol = snes->vec_sol_always = x;
938 
939   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
940   /*
941       This version replaces the user provided Jacobian matrix with a
942       matrix-free version but still employs the user-provided preconditioner matrix
943   */
944   if (flg) {
945     Mat J;
946     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
947     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
948     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n");
949     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
950     ierr = MatDestroy(J);CHKERRQ(ierr);
951   }
952 
953 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
954   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
955   if (flg) {
956     Mat J;
957     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
958     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
959     ierr = MatDestroy(J);CHKERRQ(ierr);
960   }
961 #endif
962 
963   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
964   /*
965       This version replaces both the user-provided Jacobian and the user-
966       provided preconditioner matrix with the default matrix free version.
967    */
968   if (flg) {
969     Mat  J;
970     SLES sles;
971     PC   pc;
972     KSP  ksp;
973 
974     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
975     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
976     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
977     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
978     ierr = MatDestroy(J);CHKERRQ(ierr);
979 
980     /* force no preconditioner */
981     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
982     ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
983     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
984     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
985     if (!flg) {
986       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
987     }
988   }
989 
990   if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
991   if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
992   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
993   if (snes->vec_func == snes->vec_sol) {
994     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
995   }
996 
997   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
998   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
999   if (snes->ksp_ewconv && !iseqtr) {
1000     SLES sles; KSP ksp;
1001     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1002     ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1003     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1004   }
1005 
1006   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
1007   snes->setupcalled = 1;
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "SNESDestroy"
1013 /*@C
1014    SNESDestroy - Destroys the nonlinear solver context that was created
1015    with SNESCreate().
1016 
1017    Collective on SNES
1018 
1019    Input Parameter:
1020 .  snes - the SNES context
1021 
1022    Level: beginner
1023 
1024 .keywords: SNES, nonlinear, destroy
1025 
1026 .seealso: SNESCreate(), SNESSolve()
1027 @*/
1028 int SNESDestroy(SNES snes)
1029 {
1030   int i,ierr;
1031 
1032   PetscFunctionBegin;
1033   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1034   if (--snes->refct > 0) PetscFunctionReturn(0);
1035 
1036   /* if memory was published with AMS then destroy it */
1037   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1038 
1039   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1040   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
1041   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1042   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1043   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1044   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1045   for (i=0; i<snes->numbermonitors; i++) {
1046     if (snes->monitordestroy[i]) {
1047       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1048     }
1049   }
1050   PetscLogObjectDestroy((PetscObject)snes);
1051   PetscHeaderDestroy((PetscObject)snes);
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 /* ----------- Routines to set solver parameters ---------- */
1056 
1057 #undef __FUNCT__
1058 #define __FUNCT__ "SNESSetTolerances"
1059 /*@
1060    SNESSetTolerances - Sets various parameters used in convergence tests.
1061 
1062    Collective on SNES
1063 
1064    Input Parameters:
1065 +  snes - the SNES context
1066 .  atol - absolute convergence tolerance
1067 .  rtol - relative convergence tolerance
1068 .  stol -  convergence tolerance in terms of the norm
1069            of the change in the solution between steps
1070 .  maxit - maximum number of iterations
1071 -  maxf - maximum number of function evaluations
1072 
1073    Options Database Keys:
1074 +    -snes_atol <atol> - Sets atol
1075 .    -snes_rtol <rtol> - Sets rtol
1076 .    -snes_stol <stol> - Sets stol
1077 .    -snes_max_it <maxit> - Sets maxit
1078 -    -snes_max_funcs <maxf> - Sets maxf
1079 
1080    Notes:
1081    The default maximum number of iterations is 50.
1082    The default maximum number of function evaluations is 1000.
1083 
1084    Level: intermediate
1085 
1086 .keywords: SNES, nonlinear, set, convergence, tolerances
1087 
1088 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1089 @*/
1090 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
1091 {
1092   PetscFunctionBegin;
1093   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1094   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1095   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1096   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1097   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1098   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1099   PetscFunctionReturn(0);
1100 }
1101 
1102 #undef __FUNCT__
1103 #define __FUNCT__ "SNESGetTolerances"
1104 /*@
1105    SNESGetTolerances - Gets various parameters used in convergence tests.
1106 
1107    Not Collective
1108 
1109    Input Parameters:
1110 +  snes - the SNES context
1111 .  atol - absolute convergence tolerance
1112 .  rtol - relative convergence tolerance
1113 .  stol -  convergence tolerance in terms of the norm
1114            of the change in the solution between steps
1115 .  maxit - maximum number of iterations
1116 -  maxf - maximum number of function evaluations
1117 
1118    Notes:
1119    The user can specify PETSC_NULL for any parameter that is not needed.
1120 
1121    Level: intermediate
1122 
1123 .keywords: SNES, nonlinear, get, convergence, tolerances
1124 
1125 .seealso: SNESSetTolerances()
1126 @*/
1127 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
1128 {
1129   PetscFunctionBegin;
1130   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1131   if (atol)  *atol  = snes->atol;
1132   if (rtol)  *rtol  = snes->rtol;
1133   if (stol)  *stol  = snes->xtol;
1134   if (maxit) *maxit = snes->max_its;
1135   if (maxf)  *maxf  = snes->max_funcs;
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 #undef __FUNCT__
1140 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1141 /*@
1142    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1143 
1144    Collective on SNES
1145 
1146    Input Parameters:
1147 +  snes - the SNES context
1148 -  tol - tolerance
1149 
1150    Options Database Key:
1151 .  -snes_trtol <tol> - Sets tol
1152 
1153    Level: intermediate
1154 
1155 .keywords: SNES, nonlinear, set, trust region, tolerance
1156 
1157 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1158 @*/
1159 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1160 {
1161   PetscFunctionBegin;
1162   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1163   snes->deltatol = tol;
1164   PetscFunctionReturn(0);
1165 }
1166 
1167 /*
1168    Duplicate the lg monitors for SNES from KSP; for some reason with
1169    dynamic libraries things don't work under Sun4 if we just use
1170    macros instead of functions
1171 */
1172 #undef __FUNCT__
1173 #define __FUNCT__ "SNESLGMonitor"
1174 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1175 {
1176   int ierr;
1177 
1178   PetscFunctionBegin;
1179   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1180   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 #undef __FUNCT__
1185 #define __FUNCT__ "SNESLGMonitorCreate"
1186 int SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1187 {
1188   int ierr;
1189 
1190   PetscFunctionBegin;
1191   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1192   PetscFunctionReturn(0);
1193 }
1194 
1195 #undef __FUNCT__
1196 #define __FUNCT__ "SNESLGMonitorDestroy"
1197 int SNESLGMonitorDestroy(PetscDrawLG draw)
1198 {
1199   int ierr;
1200 
1201   PetscFunctionBegin;
1202   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 /* ------------ Routines to set performance monitoring options ----------- */
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "SNESSetMonitor"
1210 /*@C
1211    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
1212    iteration of the nonlinear solver to display the iteration's
1213    progress.
1214 
1215    Collective on SNES
1216 
1217    Input Parameters:
1218 +  snes - the SNES context
1219 .  func - monitoring routine
1220 .  mctx - [optional] user-defined context for private data for the
1221           monitor routine (use PETSC_NULL if no context is desitre)
1222 -  monitordestroy - [optional] routine that frees monitor context
1223           (may be PETSC_NULL)
1224 
1225    Calling sequence of func:
1226 $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1227 
1228 +    snes - the SNES context
1229 .    its - iteration number
1230 .    norm - 2-norm function value (may be estimated)
1231 -    mctx - [optional] monitoring context
1232 
1233    Options Database Keys:
1234 +    -snes_monitor        - sets SNESDefaultMonitor()
1235 .    -snes_xmonitor       - sets line graph monitor,
1236                             uses SNESLGMonitorCreate()
1237 _    -snes_cancelmonitors - cancels all monitors that have
1238                             been hardwired into a code by
1239                             calls to SNESSetMonitor(), but
1240                             does not cancel those set via
1241                             the options database.
1242 
1243    Notes:
1244    Several different monitoring routines may be set by calling
1245    SNESSetMonitor() multiple times; all will be called in the
1246    order in which they were set.
1247 
1248    Level: intermediate
1249 
1250 .keywords: SNES, nonlinear, set, monitor
1251 
1252 .seealso: SNESDefaultMonitor(), SNESClearMonitor()
1253 @*/
1254 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
1255 {
1256   PetscFunctionBegin;
1257   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1258   if (snes->numbermonitors >= MAXSNESMONITORS) {
1259     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1260   }
1261 
1262   snes->monitor[snes->numbermonitors]           = func;
1263   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1264   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1265   PetscFunctionReturn(0);
1266 }
1267 
1268 #undef __FUNCT__
1269 #define __FUNCT__ "SNESClearMonitor"
1270 /*@C
1271    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1272 
1273    Collective on SNES
1274 
1275    Input Parameters:
1276 .  snes - the SNES context
1277 
1278    Options Database Key:
1279 .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1280     into a code by calls to SNESSetMonitor(), but does not cancel those
1281     set via the options database
1282 
1283    Notes:
1284    There is no way to clear one specific monitor from a SNES object.
1285 
1286    Level: intermediate
1287 
1288 .keywords: SNES, nonlinear, set, monitor
1289 
1290 .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1291 @*/
1292 int SNESClearMonitor(SNES snes)
1293 {
1294   PetscFunctionBegin;
1295   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1296   snes->numbermonitors = 0;
1297   PetscFunctionReturn(0);
1298 }
1299 
1300 #undef __FUNCT__
1301 #define __FUNCT__ "SNESSetConvergenceTest"
1302 /*@C
1303    SNESSetConvergenceTest - Sets the function that is to be used
1304    to test for convergence of the nonlinear iterative solution.
1305 
1306    Collective on SNES
1307 
1308    Input Parameters:
1309 +  snes - the SNES context
1310 .  func - routine to test for convergence
1311 -  cctx - [optional] context for private data for the convergence routine
1312           (may be PETSC_NULL)
1313 
1314    Calling sequence of func:
1315 $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1316 
1317 +    snes - the SNES context
1318 .    cctx - [optional] convergence context
1319 .    reason - reason for convergence/divergence
1320 .    xnorm - 2-norm of current iterate
1321 .    gnorm - 2-norm of current step
1322 -    f - 2-norm of function
1323 
1324    Level: advanced
1325 
1326 .keywords: SNES, nonlinear, set, convergence, test
1327 
1328 .seealso: SNESConverged_LS(), SNESConverged_TR()
1329 @*/
1330 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
1331 {
1332   PetscFunctionBegin;
1333   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1334   (snes)->converged = func;
1335   (snes)->cnvP      = cctx;
1336   PetscFunctionReturn(0);
1337 }
1338 
1339 #undef __FUNCT__
1340 #define __FUNCT__ "SNESGetConvergedReason"
1341 /*@C
1342    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1343 
1344    Not Collective
1345 
1346    Input Parameter:
1347 .  snes - the SNES context
1348 
1349    Output Parameter:
1350 .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1351             manual pages for the individual convergence tests for complete lists
1352 
1353    Level: intermediate
1354 
1355    Notes: Can only be called after the call the SNESSolve() is complete.
1356 
1357 .keywords: SNES, nonlinear, set, convergence, test
1358 
1359 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1360 @*/
1361 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1362 {
1363   PetscFunctionBegin;
1364   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1365   *reason = snes->reason;
1366   PetscFunctionReturn(0);
1367 }
1368 
1369 #undef __FUNCT__
1370 #define __FUNCT__ "SNESSetConvergenceHistory"
1371 /*@
1372    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1373 
1374    Collective on SNES
1375 
1376    Input Parameters:
1377 +  snes - iterative context obtained from SNESCreate()
1378 .  a   - array to hold history
1379 .  its - integer array holds the number of linear iterations for each solve.
1380 .  na  - size of a and its
1381 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1382            else it continues storing new values for new nonlinear solves after the old ones
1383 
1384    Notes:
1385    If set, this array will contain the function norms computed
1386    at each step.
1387 
1388    This routine is useful, e.g., when running a code for purposes
1389    of accurate performance monitoring, when no I/O should be done
1390    during the section of code that is being timed.
1391 
1392    Level: intermediate
1393 
1394 .keywords: SNES, set, convergence, history
1395 
1396 .seealso: SNESGetConvergenceHistory()
1397 
1398 @*/
1399 int SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset)
1400 {
1401   PetscFunctionBegin;
1402   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1403   if (na) PetscValidScalarPointer(a);
1404   snes->conv_hist       = a;
1405   snes->conv_hist_its   = its;
1406   snes->conv_hist_max   = na;
1407   snes->conv_hist_reset = reset;
1408   PetscFunctionReturn(0);
1409 }
1410 
1411 #undef __FUNCT__
1412 #define __FUNCT__ "SNESGetConvergenceHistory"
1413 /*@C
1414    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1415 
1416    Collective on SNES
1417 
1418    Input Parameter:
1419 .  snes - iterative context obtained from SNESCreate()
1420 
1421    Output Parameters:
1422 .  a   - array to hold history
1423 .  its - integer array holds the number of linear iterations (or
1424          negative if not converged) for each solve.
1425 -  na  - size of a and its
1426 
1427    Notes:
1428     The calling sequence for this routine in Fortran is
1429 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1430 
1431    This routine is useful, e.g., when running a code for purposes
1432    of accurate performance monitoring, when no I/O should be done
1433    during the section of code that is being timed.
1434 
1435    Level: intermediate
1436 
1437 .keywords: SNES, get, convergence, history
1438 
1439 .seealso: SNESSetConvergencHistory()
1440 
1441 @*/
1442 int SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na)
1443 {
1444   PetscFunctionBegin;
1445   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1446   if (a)   *a   = snes->conv_hist;
1447   if (its) *its = snes->conv_hist_its;
1448   if (na) *na   = snes->conv_hist_len;
1449   PetscFunctionReturn(0);
1450 }
1451 
1452 #undef __FUNCT__
1453 #define __FUNCT__ "SNESSetRhsBC"
1454 /*@
1455   SNESSetRhsBC - Sets the function which applies boundary conditions
1456   to the Rhs of each system.
1457 
1458   Collective on SNES
1459 
1460   Input Parameters:
1461 . snes - The nonlinear solver context
1462 . func - The function
1463 
1464   Calling sequence of func:
1465 . func (SNES snes, Vec rhs, void *ctx);
1466 
1467 . rhs - The current rhs vector
1468 . ctx - The user-context
1469 
1470   Level: intermediate
1471 
1472 .keywords: SNES, Rhs, boundary conditions
1473 .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
1474 @*/
1475 int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *))
1476 {
1477   PetscFunctionBegin;
1478   PetscValidHeaderSpecific(snes, SNES_COOKIE);
1479   snes->applyrhsbc = func;
1480   PetscFunctionReturn(0);
1481 }
1482 
1483 #undef __FUNCT__
1484 #define __FUNCT__ "SNESDefaultRhsBC"
1485 /*@
1486   SNESDefaultRhsBC - The default boundary condition function which does nothing.
1487 
1488   Not collective
1489 
1490   Input Parameters:
1491 . snes - The nonlinear solver context
1492 . rhs  - The Rhs
1493 . ctx  - The user-context
1494 
1495   Level: beginner
1496 
1497 .keywords: SNES, Rhs, boundary conditions
1498 .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
1499 @*/
1500 int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
1501 {
1502   PetscFunctionBegin;
1503   PetscFunctionReturn(0);
1504 }
1505 
1506 #undef __FUNCT__
1507 #define __FUNCT__ "SNESSetSolutionBC"
1508 /*@
1509   SNESSetSolutionBC - Sets the function which applies boundary conditions
1510   to the solution of each system.
1511 
1512   Collective on SNES
1513 
1514   Input Parameters:
1515 . snes - The nonlinear solver context
1516 . func - The function
1517 
1518   Calling sequence of func:
1519 . func (SNES snes, Vec rsol, void *ctx);
1520 
1521 . sol - The current solution vector
1522 . ctx - The user-context
1523 
1524   Level: intermediate
1525 
1526 .keywords: SNES, solution, boundary conditions
1527 .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
1528 @*/
1529 int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *))
1530 {
1531   PetscFunctionBegin;
1532   PetscValidHeaderSpecific(snes, SNES_COOKIE);
1533   snes->applysolbc = func;
1534   PetscFunctionReturn(0);
1535 }
1536 
1537 #undef __FUNCT__
1538 #define __FUNCT__ "SNESDefaultSolutionBC"
1539 /*@
1540   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
1541 
1542   Not collective
1543 
1544   Input Parameters:
1545 . snes - The nonlinear solver context
1546 . sol  - The solution
1547 . ctx  - The user-context
1548 
1549   Level: beginner
1550 
1551 .keywords: SNES, solution, boundary conditions
1552 .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
1553 @*/
1554 int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
1555 {
1556   PetscFunctionBegin;
1557   PetscFunctionReturn(0);
1558 }
1559 
1560 #undef __FUNCT__
1561 #define __FUNCT__ "SNESSetUpdate"
1562 /*@
1563   SNESSetUpdate - Sets the general-purpose update function called
1564   at the beginning of every step of the iteration.
1565 
1566   Collective on SNES
1567 
1568   Input Parameters:
1569 . snes - The nonlinear solver context
1570 . func - The function
1571 
1572   Calling sequence of func:
1573 . func (TS ts, int step);
1574 
1575 . step - The current step of the iteration
1576 
1577   Level: intermediate
1578 
1579 .keywords: SNES, update
1580 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
1581 @*/
1582 int SNESSetUpdate(SNES snes, int (*func)(SNES, int))
1583 {
1584   PetscFunctionBegin;
1585   PetscValidHeaderSpecific(snes, SNES_COOKIE);
1586   snes->update = func;
1587   PetscFunctionReturn(0);
1588 }
1589 
1590 #undef __FUNCT__
1591 #define __FUNCT__ "SNESDefaultUpdate"
1592 /*@
1593   SNESDefaultUpdate - The default update function which does nothing.
1594 
1595   Not collective
1596 
1597   Input Parameters:
1598 . snes - The nonlinear solver context
1599 . step - The current step of the iteration
1600 
1601   Level: intermediate
1602 
1603 .keywords: SNES, update
1604 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
1605 @*/
1606 int SNESDefaultUpdate(SNES snes, int step)
1607 {
1608   PetscFunctionBegin;
1609   PetscFunctionReturn(0);
1610 }
1611 
1612 #undef __FUNCT__
1613 #define __FUNCT__ "SNESScaleStep_Private"
1614 /*
1615    SNESScaleStep_Private - Scales a step so that its length is less than the
1616    positive parameter delta.
1617 
1618     Input Parameters:
1619 +   snes - the SNES context
1620 .   y - approximate solution of linear system
1621 .   fnorm - 2-norm of current function
1622 -   delta - trust region size
1623 
1624     Output Parameters:
1625 +   gpnorm - predicted function norm at the new point, assuming local
1626     linearization.  The value is zero if the step lies within the trust
1627     region, and exceeds zero otherwise.
1628 -   ynorm - 2-norm of the step
1629 
1630     Note:
1631     For non-trust region methods such as SNESLS, the parameter delta
1632     is set to be the maximum allowable step size.
1633 
1634 .keywords: SNES, nonlinear, scale, step
1635 */
1636 int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
1637 {
1638   PetscReal   nrm;
1639   PetscScalar cnorm;
1640   int         ierr;
1641 
1642   PetscFunctionBegin;
1643   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1644   PetscValidHeaderSpecific(y,VEC_COOKIE);
1645   PetscCheckSameComm(snes,y);
1646 
1647   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1648   if (nrm > *delta) {
1649      nrm = *delta/nrm;
1650      *gpnorm = (1.0 - nrm)*(*fnorm);
1651      cnorm = nrm;
1652      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
1653      *ynorm = *delta;
1654   } else {
1655      *gpnorm = 0.0;
1656      *ynorm = nrm;
1657   }
1658   PetscFunctionReturn(0);
1659 }
1660 
1661 #undef __FUNCT__
1662 #define __FUNCT__ "SNESSolve"
1663 /*@
1664    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1665    SNESCreate() and optional routines of the form SNESSetXXX().
1666 
1667    Collective on SNES
1668 
1669    Input Parameters:
1670 +  snes - the SNES context
1671 -  x - the solution vector
1672 
1673    Output Parameter:
1674 .  its - number of iterations until termination
1675 
1676    Notes:
1677    The user should initialize the vector,x, with the initial guess
1678    for the nonlinear solve prior to calling SNESSolve.  In particular,
1679    to employ an initial guess of zero, the user should explicitly set
1680    this vector to zero by calling VecSet().
1681 
1682    Level: beginner
1683 
1684 .keywords: SNES, nonlinear, solve
1685 
1686 .seealso: SNESCreate(), SNESDestroy()
1687 @*/
1688 int SNESSolve(SNES snes,Vec x,int *its)
1689 {
1690   int        ierr;
1691   PetscTruth flg;
1692 
1693   PetscFunctionBegin;
1694   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1695   PetscValidHeaderSpecific(x,VEC_COOKIE);
1696   PetscCheckSameComm(snes,x);
1697   PetscValidIntPointer(its);
1698   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
1699 
1700   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
1701   else {snes->vec_sol = snes->vec_sol_always = x;}
1702   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
1703   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1704   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1705   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
1706   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1707   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
1708   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); }
1709   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1710   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
1711   PetscFunctionReturn(0);
1712 }
1713 
1714 /* --------- Internal routines for SNES Package --------- */
1715 
1716 #undef __FUNCT__
1717 #define __FUNCT__ "SNESSetType"
1718 /*@C
1719    SNESSetType - Sets the method for the nonlinear solver.
1720 
1721    Collective on SNES
1722 
1723    Input Parameters:
1724 +  snes - the SNES context
1725 -  type - a known method
1726 
1727    Options Database Key:
1728 .  -snes_type <type> - Sets the method; use -help for a list
1729    of available methods (for instance, ls or tr)
1730 
1731    Notes:
1732    See "petsc/include/petscsnes.h" for available methods (for instance)
1733 +    SNESLS - Newton's method with line search
1734      (systems of nonlinear equations)
1735 .    SNESTR - Newton's method with trust region
1736      (systems of nonlinear equations)
1737 
1738   Normally, it is best to use the SNESSetFromOptions() command and then
1739   set the SNES solver type from the options database rather than by using
1740   this routine.  Using the options database provides the user with
1741   maximum flexibility in evaluating the many nonlinear solvers.
1742   The SNESSetType() routine is provided for those situations where it
1743   is necessary to set the nonlinear solver independently of the command
1744   line or options database.  This might be the case, for example, when
1745   the choice of solver changes during the execution of the program,
1746   and the user's application is taking responsibility for choosing the
1747   appropriate method.
1748 
1749   Level: intermediate
1750 
1751 .keywords: SNES, set, type
1752 
1753 .seealso: SNESType, SNESCreate()
1754 
1755 @*/
1756 int SNESSetType(SNES snes,SNESType type)
1757 {
1758   int        ierr,(*r)(SNES);
1759   PetscTruth match;
1760 
1761   PetscFunctionBegin;
1762   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1763   PetscValidCharPointer(type);
1764 
1765   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
1766   if (match) PetscFunctionReturn(0);
1767 
1768   if (snes->setupcalled) {
1769     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
1770     snes->data = 0;
1771   }
1772 
1773   /* Get the function pointers for the iterative method requested */
1774   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1775 
1776   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1777 
1778   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
1779 
1780   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
1781   snes->data = 0;
1782   ierr = (*r)(snes);CHKERRQ(ierr);
1783   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
1784 
1785   PetscFunctionReturn(0);
1786 }
1787 
1788 
1789 /* --------------------------------------------------------------------- */
1790 #undef __FUNCT__
1791 #define __FUNCT__ "SNESRegisterDestroy"
1792 /*@C
1793    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1794    registered by SNESRegisterDynamic().
1795 
1796    Not Collective
1797 
1798    Level: advanced
1799 
1800 .keywords: SNES, nonlinear, register, destroy
1801 
1802 .seealso: SNESRegisterAll(), SNESRegisterAll()
1803 @*/
1804 int SNESRegisterDestroy(void)
1805 {
1806   int ierr;
1807 
1808   PetscFunctionBegin;
1809   if (SNESList) {
1810     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
1811     SNESList = 0;
1812   }
1813   SNESRegisterAllCalled = PETSC_FALSE;
1814   PetscFunctionReturn(0);
1815 }
1816 
1817 #undef __FUNCT__
1818 #define __FUNCT__ "SNESGetType"
1819 /*@C
1820    SNESGetType - Gets the SNES method type and name (as a string).
1821 
1822    Not Collective
1823 
1824    Input Parameter:
1825 .  snes - nonlinear solver context
1826 
1827    Output Parameter:
1828 .  type - SNES method (a character string)
1829 
1830    Level: intermediate
1831 
1832 .keywords: SNES, nonlinear, get, type, name
1833 @*/
1834 int SNESGetType(SNES snes,SNESType *type)
1835 {
1836   PetscFunctionBegin;
1837   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1838   *type = snes->type_name;
1839   PetscFunctionReturn(0);
1840 }
1841 
1842 #undef __FUNCT__
1843 #define __FUNCT__ "SNESGetSolution"
1844 /*@C
1845    SNESGetSolution - Returns the vector where the approximate solution is
1846    stored.
1847 
1848    Not Collective, but Vec is parallel if SNES is parallel
1849 
1850    Input Parameter:
1851 .  snes - the SNES context
1852 
1853    Output Parameter:
1854 .  x - the solution
1855 
1856    Level: advanced
1857 
1858 .keywords: SNES, nonlinear, get, solution
1859 
1860 .seealso: SNESGetFunction(), SNESGetSolutionUpdate()
1861 @*/
1862 int SNESGetSolution(SNES snes,Vec *x)
1863 {
1864   PetscFunctionBegin;
1865   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1866   *x = snes->vec_sol_always;
1867   PetscFunctionReturn(0);
1868 }
1869 
1870 #undef __FUNCT__
1871 #define __FUNCT__ "SNESGetSolutionUpdate"
1872 /*@C
1873    SNESGetSolutionUpdate - Returns the vector where the solution update is
1874    stored.
1875 
1876    Not Collective, but Vec is parallel if SNES is parallel
1877 
1878    Input Parameter:
1879 .  snes - the SNES context
1880 
1881    Output Parameter:
1882 .  x - the solution update
1883 
1884    Level: advanced
1885 
1886 .keywords: SNES, nonlinear, get, solution, update
1887 
1888 .seealso: SNESGetSolution(), SNESGetFunction
1889 @*/
1890 int SNESGetSolutionUpdate(SNES snes,Vec *x)
1891 {
1892   PetscFunctionBegin;
1893   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1894   *x = snes->vec_sol_update_always;
1895   PetscFunctionReturn(0);
1896 }
1897 
1898 #undef __FUNCT__
1899 #define __FUNCT__ "SNESGetFunction"
1900 /*@C
1901    SNESGetFunction - Returns the vector where the function is stored.
1902 
1903    Not Collective, but Vec is parallel if SNES is parallel
1904 
1905    Input Parameter:
1906 .  snes - the SNES context
1907 
1908    Output Parameter:
1909 +  r - the function (or PETSC_NULL)
1910 .  ctx - the function context (or PETSC_NULL)
1911 -  func - the function (or PETSC_NULL)
1912 
1913    Level: advanced
1914 
1915 .keywords: SNES, nonlinear, get, function
1916 
1917 .seealso: SNESSetFunction(), SNESGetSolution()
1918 @*/
1919 int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
1920 {
1921   PetscFunctionBegin;
1922   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1923   if (r)    *r    = snes->vec_func_always;
1924   if (ctx)  *ctx  = snes->funP;
1925   if (func) *func = snes->computefunction;
1926   PetscFunctionReturn(0);
1927 }
1928 
1929 #undef __FUNCT__
1930 #define __FUNCT__ "SNESSetOptionsPrefix"
1931 /*@C
1932    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1933    SNES options in the database.
1934 
1935    Collective on SNES
1936 
1937    Input Parameter:
1938 +  snes - the SNES context
1939 -  prefix - the prefix to prepend to all option names
1940 
1941    Notes:
1942    A hyphen (-) must NOT be given at the beginning of the prefix name.
1943    The first character of all runtime options is AUTOMATICALLY the hyphen.
1944 
1945    Level: advanced
1946 
1947 .keywords: SNES, set, options, prefix, database
1948 
1949 .seealso: SNESSetFromOptions()
1950 @*/
1951 int SNESSetOptionsPrefix(SNES snes,const char prefix[])
1952 {
1953   int ierr;
1954 
1955   PetscFunctionBegin;
1956   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1957   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
1958   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1959   PetscFunctionReturn(0);
1960 }
1961 
1962 #undef __FUNCT__
1963 #define __FUNCT__ "SNESAppendOptionsPrefix"
1964 /*@C
1965    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1966    SNES options in the database.
1967 
1968    Collective on SNES
1969 
1970    Input Parameters:
1971 +  snes - the SNES context
1972 -  prefix - the prefix to prepend to all option names
1973 
1974    Notes:
1975    A hyphen (-) must NOT be given at the beginning of the prefix name.
1976    The first character of all runtime options is AUTOMATICALLY the hyphen.
1977 
1978    Level: advanced
1979 
1980 .keywords: SNES, append, options, prefix, database
1981 
1982 .seealso: SNESGetOptionsPrefix()
1983 @*/
1984 int SNESAppendOptionsPrefix(SNES snes,const char prefix[])
1985 {
1986   int ierr;
1987 
1988   PetscFunctionBegin;
1989   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1990   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
1991   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1992   PetscFunctionReturn(0);
1993 }
1994 
1995 #undef __FUNCT__
1996 #define __FUNCT__ "SNESGetOptionsPrefix"
1997 /*@C
1998    SNESGetOptionsPrefix - Sets the prefix used for searching for all
1999    SNES options in the database.
2000 
2001    Not Collective
2002 
2003    Input Parameter:
2004 .  snes - the SNES context
2005 
2006    Output Parameter:
2007 .  prefix - pointer to the prefix string used
2008 
2009    Notes: On the fortran side, the user should pass in a string 'prifix' of
2010    sufficient length to hold the prefix.
2011 
2012    Level: advanced
2013 
2014 .keywords: SNES, get, options, prefix, database
2015 
2016 .seealso: SNESAppendOptionsPrefix()
2017 @*/
2018 int SNESGetOptionsPrefix(SNES snes,char *prefix[])
2019 {
2020   int ierr;
2021 
2022   PetscFunctionBegin;
2023   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2024   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2025   PetscFunctionReturn(0);
2026 }
2027 
2028 
2029 #undef __FUNCT__
2030 #define __FUNCT__ "SNESRegister"
2031 /*@C
2032   SNESRegister - See SNESRegisterDynamic()
2033 
2034   Level: advanced
2035 @*/
2036 int SNESRegister(const char sname[],const char path[],const char name[],int (*function)(SNES))
2037 {
2038   char fullname[256];
2039   int  ierr;
2040 
2041   PetscFunctionBegin;
2042   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2043   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2044   PetscFunctionReturn(0);
2045 }
2046 
2047 #undef __FUNCT__
2048 #define __FUNCT__ "SNESTestLocalMin"
2049 int SNESTestLocalMin(SNES snes)
2050 {
2051   int         ierr,N,i,j;
2052   Vec         u,uh,fh;
2053   PetscScalar value;
2054   PetscReal   norm;
2055 
2056   PetscFunctionBegin;
2057   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2058   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2059   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2060 
2061   /* currently only works for sequential */
2062   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2063   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2064   for (i=0; i<N; i++) {
2065     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2066     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr);
2067     for (j=-10; j<11; j++) {
2068       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2069       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2070       ierr  = (*snes->computefunction)(snes,uh,fh,snes->funP);CHKERRQ(ierr);
2071       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2072       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %d %18.16e\n",j,norm);CHKERRQ(ierr);
2073       value = -value;
2074       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2075     }
2076   }
2077   ierr = VecDestroy(uh);CHKERRQ(ierr);
2078   ierr = VecDestroy(fh);CHKERRQ(ierr);
2079   PetscFunctionReturn(0);
2080 }
2081