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