xref: /petsc/src/snes/interface/snes.c (revision 2c62819748cfd70e832b5b508eca90f58e6249c6)
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 SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval;
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: SLESGetPC(), 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 
973     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
974     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
975     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
976     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
977     ierr = MatDestroy(J);CHKERRQ(ierr);
978 
979     /* force no preconditioner */
980     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
981     ierr = SLESGetPC(sles,&pc);CHKERRQ(ierr);
982     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
983   }
984 
985   if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
986   if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
987   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
988   if (snes->vec_func == snes->vec_sol) {
989     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
990   }
991 
992   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
993   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
994   if (snes->ksp_ewconv && !iseqtr) {
995     SLES sles; KSP ksp;
996     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
997     ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
998     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
999   }
1000 
1001   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
1002   snes->setupcalled = 1;
1003   PetscFunctionReturn(0);
1004 }
1005 
1006 #undef __FUNCT__
1007 #define __FUNCT__ "SNESDestroy"
1008 /*@C
1009    SNESDestroy - Destroys the nonlinear solver context that was created
1010    with SNESCreate().
1011 
1012    Collective on SNES
1013 
1014    Input Parameter:
1015 .  snes - the SNES context
1016 
1017    Level: beginner
1018 
1019 .keywords: SNES, nonlinear, destroy
1020 
1021 .seealso: SNESCreate(), SNESSolve()
1022 @*/
1023 int SNESDestroy(SNES snes)
1024 {
1025   int i,ierr;
1026 
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1029   if (--snes->refct > 0) PetscFunctionReturn(0);
1030 
1031   /* if memory was published with AMS then destroy it */
1032   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1033 
1034   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1035   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
1036   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1037   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1038   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1039   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1040   for (i=0; i<snes->numbermonitors; i++) {
1041     if (snes->monitordestroy[i]) {
1042       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1043     }
1044   }
1045   PetscLogObjectDestroy((PetscObject)snes);
1046   PetscHeaderDestroy((PetscObject)snes);
1047   PetscFunctionReturn(0);
1048 }
1049 
1050 /* ----------- Routines to set solver parameters ---------- */
1051 
1052 #undef __FUNCT__
1053 #define __FUNCT__ "SNESSetTolerances"
1054 /*@
1055    SNESSetTolerances - Sets various parameters used in convergence tests.
1056 
1057    Collective on SNES
1058 
1059    Input Parameters:
1060 +  snes - the SNES context
1061 .  atol - absolute convergence tolerance
1062 .  rtol - relative convergence tolerance
1063 .  stol -  convergence tolerance in terms of the norm
1064            of the change in the solution between steps
1065 .  maxit - maximum number of iterations
1066 -  maxf - maximum number of function evaluations
1067 
1068    Options Database Keys:
1069 +    -snes_atol <atol> - Sets atol
1070 .    -snes_rtol <rtol> - Sets rtol
1071 .    -snes_stol <stol> - Sets stol
1072 .    -snes_max_it <maxit> - Sets maxit
1073 -    -snes_max_funcs <maxf> - Sets maxf
1074 
1075    Notes:
1076    The default maximum number of iterations is 50.
1077    The default maximum number of function evaluations is 1000.
1078 
1079    Level: intermediate
1080 
1081 .keywords: SNES, nonlinear, set, convergence, tolerances
1082 
1083 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1084 @*/
1085 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
1086 {
1087   PetscFunctionBegin;
1088   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1089   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1090   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1091   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1092   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1093   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1094   PetscFunctionReturn(0);
1095 }
1096 
1097 #undef __FUNCT__
1098 #define __FUNCT__ "SNESGetTolerances"
1099 /*@
1100    SNESGetTolerances - Gets various parameters used in convergence tests.
1101 
1102    Not Collective
1103 
1104    Input Parameters:
1105 +  snes - the SNES context
1106 .  atol - absolute convergence tolerance
1107 .  rtol - relative convergence tolerance
1108 .  stol -  convergence tolerance in terms of the norm
1109            of the change in the solution between steps
1110 .  maxit - maximum number of iterations
1111 -  maxf - maximum number of function evaluations
1112 
1113    Notes:
1114    The user can specify PETSC_NULL for any parameter that is not needed.
1115 
1116    Level: intermediate
1117 
1118 .keywords: SNES, nonlinear, get, convergence, tolerances
1119 
1120 .seealso: SNESSetTolerances()
1121 @*/
1122 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
1123 {
1124   PetscFunctionBegin;
1125   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1126   if (atol)  *atol  = snes->atol;
1127   if (rtol)  *rtol  = snes->rtol;
1128   if (stol)  *stol  = snes->xtol;
1129   if (maxit) *maxit = snes->max_its;
1130   if (maxf)  *maxf  = snes->max_funcs;
1131   PetscFunctionReturn(0);
1132 }
1133 
1134 #undef __FUNCT__
1135 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1136 /*@
1137    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1138 
1139    Collective on SNES
1140 
1141    Input Parameters:
1142 +  snes - the SNES context
1143 -  tol - tolerance
1144 
1145    Options Database Key:
1146 .  -snes_trtol <tol> - Sets tol
1147 
1148    Level: intermediate
1149 
1150 .keywords: SNES, nonlinear, set, trust region, tolerance
1151 
1152 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1153 @*/
1154 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1155 {
1156   PetscFunctionBegin;
1157   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1158   snes->deltatol = tol;
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 /*
1163    Duplicate the lg monitors for SNES from KSP; for some reason with
1164    dynamic libraries things don't work under Sun4 if we just use
1165    macros instead of functions
1166 */
1167 #undef __FUNCT__
1168 #define __FUNCT__ "SNESLGMonitor"
1169 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1170 {
1171   int ierr;
1172 
1173   PetscFunctionBegin;
1174   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1175   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1176   PetscFunctionReturn(0);
1177 }
1178 
1179 #undef __FUNCT__
1180 #define __FUNCT__ "SNESLGMonitorCreate"
1181 int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1182 {
1183   int ierr;
1184 
1185   PetscFunctionBegin;
1186   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1187   PetscFunctionReturn(0);
1188 }
1189 
1190 #undef __FUNCT__
1191 #define __FUNCT__ "SNESLGMonitorDestroy"
1192 int SNESLGMonitorDestroy(PetscDrawLG draw)
1193 {
1194   int ierr;
1195 
1196   PetscFunctionBegin;
1197   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1198   PetscFunctionReturn(0);
1199 }
1200 
1201 /* ------------ Routines to set performance monitoring options ----------- */
1202 
1203 #undef __FUNCT__
1204 #define __FUNCT__ "SNESSetMonitor"
1205 /*@C
1206    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
1207    iteration of the nonlinear solver to display the iteration's
1208    progress.
1209 
1210    Collective on SNES
1211 
1212    Input Parameters:
1213 +  snes - the SNES context
1214 .  func - monitoring routine
1215 .  mctx - [optional] user-defined context for private data for the
1216           monitor routine (use PETSC_NULL if no context is desitre)
1217 -  monitordestroy - [optional] routine that frees monitor context
1218           (may be PETSC_NULL)
1219 
1220    Calling sequence of func:
1221 $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1222 
1223 +    snes - the SNES context
1224 .    its - iteration number
1225 .    norm - 2-norm function value (may be estimated)
1226 -    mctx - [optional] monitoring context
1227 
1228    Options Database Keys:
1229 +    -snes_monitor        - sets SNESDefaultMonitor()
1230 .    -snes_xmonitor       - sets line graph monitor,
1231                             uses SNESLGMonitorCreate()
1232 _    -snes_cancelmonitors - cancels all monitors that have
1233                             been hardwired into a code by
1234                             calls to SNESSetMonitor(), but
1235                             does not cancel those set via
1236                             the options database.
1237 
1238    Notes:
1239    Several different monitoring routines may be set by calling
1240    SNESSetMonitor() multiple times; all will be called in the
1241    order in which they were set.
1242 
1243    Level: intermediate
1244 
1245 .keywords: SNES, nonlinear, set, monitor
1246 
1247 .seealso: SNESDefaultMonitor(), SNESClearMonitor()
1248 @*/
1249 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
1250 {
1251   PetscFunctionBegin;
1252   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1253   if (snes->numbermonitors >= MAXSNESMONITORS) {
1254     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1255   }
1256 
1257   snes->monitor[snes->numbermonitors]           = func;
1258   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1259   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1260   PetscFunctionReturn(0);
1261 }
1262 
1263 #undef __FUNCT__
1264 #define __FUNCT__ "SNESClearMonitor"
1265 /*@C
1266    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1267 
1268    Collective on SNES
1269 
1270    Input Parameters:
1271 .  snes - the SNES context
1272 
1273    Options Database:
1274 .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1275     into a code by calls to SNESSetMonitor(), but does not cancel those
1276     set via the options database
1277 
1278    Notes:
1279    There is no way to clear one specific monitor from a SNES object.
1280 
1281    Level: intermediate
1282 
1283 .keywords: SNES, nonlinear, set, monitor
1284 
1285 .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1286 @*/
1287 int SNESClearMonitor(SNES snes)
1288 {
1289   PetscFunctionBegin;
1290   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1291   snes->numbermonitors = 0;
1292   PetscFunctionReturn(0);
1293 }
1294 
1295 #undef __FUNCT__
1296 #define __FUNCT__ "SNESSetConvergenceTest"
1297 /*@C
1298    SNESSetConvergenceTest - Sets the function that is to be used
1299    to test for convergence of the nonlinear iterative solution.
1300 
1301    Collective on SNES
1302 
1303    Input Parameters:
1304 +  snes - the SNES context
1305 .  func - routine to test for convergence
1306 -  cctx - [optional] context for private data for the convergence routine
1307           (may be PETSC_NULL)
1308 
1309    Calling sequence of func:
1310 $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1311 
1312 +    snes - the SNES context
1313 .    cctx - [optional] convergence context
1314 .    reason - reason for convergence/divergence
1315 .    xnorm - 2-norm of current iterate
1316 .    gnorm - 2-norm of current step
1317 -    f - 2-norm of function
1318 
1319    Level: advanced
1320 
1321 .keywords: SNES, nonlinear, set, convergence, test
1322 
1323 .seealso: SNESConverged_LS(), SNESConverged_TR()
1324 @*/
1325 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
1326 {
1327   PetscFunctionBegin;
1328   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1329   (snes)->converged = func;
1330   (snes)->cnvP      = cctx;
1331   PetscFunctionReturn(0);
1332 }
1333 
1334 #undef __FUNCT__
1335 #define __FUNCT__ "SNESGetConvergedReason"
1336 /*@C
1337    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1338 
1339    Not Collective
1340 
1341    Input Parameter:
1342 .  snes - the SNES context
1343 
1344    Output Parameter:
1345 .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1346             manual pages for the individual convergence tests for complete lists
1347 
1348    Level: intermediate
1349 
1350    Notes: Can only be called after the call the SNESSolve() is complete.
1351 
1352 .keywords: SNES, nonlinear, set, convergence, test
1353 
1354 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1355 @*/
1356 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1357 {
1358   PetscFunctionBegin;
1359   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1360   *reason = snes->reason;
1361   PetscFunctionReturn(0);
1362 }
1363 
1364 #undef __FUNCT__
1365 #define __FUNCT__ "SNESSetConvergenceHistory"
1366 /*@
1367    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1368 
1369    Collective on SNES
1370 
1371    Input Parameters:
1372 +  snes - iterative context obtained from SNESCreate()
1373 .  a   - array to hold history
1374 .  its - integer array holds the number of linear iterations for each solve.
1375 .  na  - size of a and its
1376 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1377            else it continues storing new values for new nonlinear solves after the old ones
1378 
1379    Notes:
1380    If set, this array will contain the function norms computed
1381    at each step.
1382 
1383    This routine is useful, e.g., when running a code for purposes
1384    of accurate performance monitoring, when no I/O should be done
1385    during the section of code that is being timed.
1386 
1387    Level: intermediate
1388 
1389 .keywords: SNES, set, convergence, history
1390 
1391 .seealso: SNESGetConvergenceHistory()
1392 
1393 @*/
1394 int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset)
1395 {
1396   PetscFunctionBegin;
1397   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1398   if (na) PetscValidScalarPointer(a);
1399   snes->conv_hist       = a;
1400   snes->conv_hist_its   = its;
1401   snes->conv_hist_max   = na;
1402   snes->conv_hist_reset = reset;
1403   PetscFunctionReturn(0);
1404 }
1405 
1406 #undef __FUNCT__
1407 #define __FUNCT__ "SNESGetConvergenceHistory"
1408 /*@C
1409    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1410 
1411    Collective on SNES
1412 
1413    Input Parameter:
1414 .  snes - iterative context obtained from SNESCreate()
1415 
1416    Output Parameters:
1417 .  a   - array to hold history
1418 .  its - integer array holds the number of linear iterations (or
1419          negative if not converged) for each solve.
1420 -  na  - size of a and its
1421 
1422    Notes:
1423     The calling sequence for this routine in Fortran is
1424 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1425 
1426    This routine is useful, e.g., when running a code for purposes
1427    of accurate performance monitoring, when no I/O should be done
1428    during the section of code that is being timed.
1429 
1430    Level: intermediate
1431 
1432 .keywords: SNES, get, convergence, history
1433 
1434 .seealso: SNESSetConvergencHistory()
1435 
1436 @*/
1437 int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na)
1438 {
1439   PetscFunctionBegin;
1440   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1441   if (a)   *a   = snes->conv_hist;
1442   if (its) *its = snes->conv_hist_its;
1443   if (na) *na   = snes->conv_hist_len;
1444   PetscFunctionReturn(0);
1445 }
1446 
1447 #undef __FUNCT__
1448 #define __FUNCT__ "SNESSetRhsBC"
1449 /*@
1450   SNESSetRhsBC - Sets the function which applies boundary conditions
1451   to the Rhs of each system.
1452 
1453   Collective on SNES
1454 
1455   Input Parameters:
1456 . snes - The nonlinear solver context
1457 . func - The function
1458 
1459   Calling sequence of func:
1460 . func (SNES snes, Vec rhs, void *ctx);
1461 
1462 . rhs - The current rhs vector
1463 . ctx - The user-context
1464 
1465   Level: intermediate
1466 
1467 .keywords: SNES, Rhs, boundary conditions
1468 .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
1469 @*/
1470 int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *))
1471 {
1472   PetscFunctionBegin;
1473   PetscValidHeaderSpecific(snes, SNES_COOKIE);
1474   snes->applyrhsbc = func;
1475   PetscFunctionReturn(0);
1476 }
1477 
1478 #undef __FUNCT__
1479 #define __FUNCT__ "SNESDefaultRhsBC"
1480 /*@
1481   SNESDefaultRhsBC - The default boundary condition function which does nothing.
1482 
1483   Not collective
1484 
1485   Input Parameters:
1486 . snes - The nonlinear solver context
1487 . rhs  - The Rhs
1488 . ctx  - The user-context
1489 
1490   Level: beginner
1491 
1492 .keywords: SNES, Rhs, boundary conditions
1493 .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
1494 @*/
1495 int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
1496 {
1497   PetscFunctionBegin;
1498   PetscFunctionReturn(0);
1499 }
1500 
1501 #undef __FUNCT__
1502 #define __FUNCT__ "SNESSetSolutionBC"
1503 /*@
1504   SNESSetSolutionBC - Sets the function which applies boundary conditions
1505   to the solution of each system.
1506 
1507   Collective on SNES
1508 
1509   Input Parameters:
1510 . snes - The nonlinear solver context
1511 . func - The function
1512 
1513   Calling sequence of func:
1514 . func (TS ts, Vec rsol, void *ctx);
1515 
1516 . sol - The current solution vector
1517 . ctx - The user-context
1518 
1519   Level: intermediate
1520 
1521 .keywords: SNES, solution, boundary conditions
1522 .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
1523 @*/
1524 int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *))
1525 {
1526   PetscFunctionBegin;
1527   PetscValidHeaderSpecific(snes, SNES_COOKIE);
1528   snes->applysolbc = func;
1529   PetscFunctionReturn(0);
1530 }
1531 
1532 #undef __FUNCT__
1533 #define __FUNCT__ "SNESDefaultSolutionBC"
1534 /*@
1535   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
1536 
1537   Not collective
1538 
1539   Input Parameters:
1540 . snes - The nonlinear solver context
1541 . sol  - The solution
1542 . ctx  - The user-context
1543 
1544   Level: beginner
1545 
1546 .keywords: SNES, solution, boundary conditions
1547 .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
1548 @*/
1549 int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
1550 {
1551   PetscFunctionBegin;
1552   PetscFunctionReturn(0);
1553 }
1554 
1555 #undef __FUNCT__
1556 #define __FUNCT__ "SNESSetUpdate"
1557 /*@
1558   SNESSetUpdate - Sets the general-purpose update function called
1559   at the beginning of every step of the iteration.
1560 
1561   Collective on SNES
1562 
1563   Input Parameters:
1564 . snes - The nonlinear solver context
1565 . func - The function
1566 
1567   Calling sequence of func:
1568 . func (TS ts, int step);
1569 
1570 . step - The current step of the iteration
1571 
1572   Level: intermediate
1573 
1574 .keywords: SNES, update
1575 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
1576 @*/
1577 int SNESSetUpdate(SNES snes, int (*func)(SNES, int))
1578 {
1579   PetscFunctionBegin;
1580   PetscValidHeaderSpecific(snes, SNES_COOKIE);
1581   snes->update = func;
1582   PetscFunctionReturn(0);
1583 }
1584 
1585 #undef __FUNCT__
1586 #define __FUNCT__ "SNESDefaultUpdate"
1587 /*@
1588   SNESDefaultUpdate - The default update function which does nothing.
1589 
1590   Not collective
1591 
1592   Input Parameters:
1593 . snes - The nonlinear solver context
1594 . step - The current step of the iteration
1595 
1596   Level: intermediate
1597 
1598 .keywords: SNES, update
1599 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
1600 @*/
1601 int SNESDefaultUpdate(SNES snes, int step)
1602 {
1603   PetscFunctionBegin;
1604   PetscFunctionReturn(0);
1605 }
1606 
1607 #undef __FUNCT__
1608 #define __FUNCT__ "SNESScaleStep_Private"
1609 /*
1610    SNESScaleStep_Private - Scales a step so that its length is less than the
1611    positive parameter delta.
1612 
1613     Input Parameters:
1614 +   snes - the SNES context
1615 .   y - approximate solution of linear system
1616 .   fnorm - 2-norm of current function
1617 -   delta - trust region size
1618 
1619     Output Parameters:
1620 +   gpnorm - predicted function norm at the new point, assuming local
1621     linearization.  The value is zero if the step lies within the trust
1622     region, and exceeds zero otherwise.
1623 -   ynorm - 2-norm of the step
1624 
1625     Note:
1626     For non-trust region methods such as SNESLS, the parameter delta
1627     is set to be the maximum allowable step size.
1628 
1629 .keywords: SNES, nonlinear, scale, step
1630 */
1631 int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
1632 {
1633   PetscReal   nrm;
1634   PetscScalar cnorm;
1635   int         ierr;
1636 
1637   PetscFunctionBegin;
1638   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1639   PetscValidHeaderSpecific(y,VEC_COOKIE);
1640   PetscCheckSameComm(snes,y);
1641 
1642   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1643   if (nrm > *delta) {
1644      nrm = *delta/nrm;
1645      *gpnorm = (1.0 - nrm)*(*fnorm);
1646      cnorm = nrm;
1647      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
1648      *ynorm = *delta;
1649   } else {
1650      *gpnorm = 0.0;
1651      *ynorm = nrm;
1652   }
1653   PetscFunctionReturn(0);
1654 }
1655 
1656 #undef __FUNCT__
1657 #define __FUNCT__ "SNESSolve"
1658 /*@
1659    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1660    SNESCreate() and optional routines of the form SNESSetXXX().
1661 
1662    Collective on SNES
1663 
1664    Input Parameters:
1665 +  snes - the SNES context
1666 -  x - the solution vector
1667 
1668    Output Parameter:
1669 .  its - number of iterations until termination
1670 
1671    Notes:
1672    The user should initialize the vector,x, with the initial guess
1673    for the nonlinear solve prior to calling SNESSolve.  In particular,
1674    to employ an initial guess of zero, the user should explicitly set
1675    this vector to zero by calling VecSet().
1676 
1677    Level: beginner
1678 
1679 .keywords: SNES, nonlinear, solve
1680 
1681 .seealso: SNESCreate(), SNESDestroy()
1682 @*/
1683 int SNESSolve(SNES snes,Vec x,int *its)
1684 {
1685   int        ierr;
1686   PetscTruth flg;
1687 
1688   PetscFunctionBegin;
1689   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1690   PetscValidHeaderSpecific(x,VEC_COOKIE);
1691   PetscCheckSameComm(snes,x);
1692   PetscValidIntPointer(its);
1693   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
1694 
1695   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
1696   else {snes->vec_sol = snes->vec_sol_always = x;}
1697   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
1698   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1699   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1700   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
1701   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1702   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
1703   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); }
1704   PetscFunctionReturn(0);
1705 }
1706 
1707 /* --------- Internal routines for SNES Package --------- */
1708 
1709 #undef __FUNCT__
1710 #define __FUNCT__ "SNESSetType"
1711 /*@C
1712    SNESSetType - Sets the method for the nonlinear solver.
1713 
1714    Collective on SNES
1715 
1716    Input Parameters:
1717 +  snes - the SNES context
1718 -  type - a known method
1719 
1720    Options Database Key:
1721 .  -snes_type <type> - Sets the method; use -help for a list
1722    of available methods (for instance, ls or tr)
1723 
1724    Notes:
1725    See "petsc/include/petscsnes.h" for available methods (for instance)
1726 +    SNESLS - Newton's method with line search
1727      (systems of nonlinear equations)
1728 .    SNESTR - Newton's method with trust region
1729      (systems of nonlinear equations)
1730 
1731   Normally, it is best to use the SNESSetFromOptions() command and then
1732   set the SNES solver type from the options database rather than by using
1733   this routine.  Using the options database provides the user with
1734   maximum flexibility in evaluating the many nonlinear solvers.
1735   The SNESSetType() routine is provided for those situations where it
1736   is necessary to set the nonlinear solver independently of the command
1737   line or options database.  This might be the case, for example, when
1738   the choice of solver changes during the execution of the program,
1739   and the user's application is taking responsibility for choosing the
1740   appropriate method.
1741 
1742   Level: intermediate
1743 
1744 .keywords: SNES, set, type
1745 
1746 .seealso: SNESType, SNESCreate()
1747 
1748 @*/
1749 int SNESSetType(SNES snes,SNESType type)
1750 {
1751   int        ierr,(*r)(SNES);
1752   PetscTruth match;
1753 
1754   PetscFunctionBegin;
1755   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1756   PetscValidCharPointer(type);
1757 
1758   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
1759   if (match) PetscFunctionReturn(0);
1760 
1761   if (snes->setupcalled) {
1762     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
1763     snes->data = 0;
1764   }
1765 
1766   /* Get the function pointers for the iterative method requested */
1767   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1768 
1769   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1770 
1771   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
1772 
1773   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
1774   snes->data = 0;
1775   ierr = (*r)(snes);CHKERRQ(ierr);
1776   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
1777 
1778   PetscFunctionReturn(0);
1779 }
1780 
1781 
1782 /* --------------------------------------------------------------------- */
1783 #undef __FUNCT__
1784 #define __FUNCT__ "SNESRegisterDestroy"
1785 /*@C
1786    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1787    registered by SNESRegisterDynamic().
1788 
1789    Not Collective
1790 
1791    Level: advanced
1792 
1793 .keywords: SNES, nonlinear, register, destroy
1794 
1795 .seealso: SNESRegisterAll(), SNESRegisterAll()
1796 @*/
1797 int SNESRegisterDestroy(void)
1798 {
1799   int ierr;
1800 
1801   PetscFunctionBegin;
1802   if (SNESList) {
1803     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
1804     SNESList = 0;
1805   }
1806   SNESRegisterAllCalled = PETSC_FALSE;
1807   PetscFunctionReturn(0);
1808 }
1809 
1810 #undef __FUNCT__
1811 #define __FUNCT__ "SNESGetType"
1812 /*@C
1813    SNESGetType - Gets the SNES method type and name (as a string).
1814 
1815    Not Collective
1816 
1817    Input Parameter:
1818 .  snes - nonlinear solver context
1819 
1820    Output Parameter:
1821 .  type - SNES method (a character string)
1822 
1823    Level: intermediate
1824 
1825 .keywords: SNES, nonlinear, get, type, name
1826 @*/
1827 int SNESGetType(SNES snes,SNESType *type)
1828 {
1829   PetscFunctionBegin;
1830   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1831   *type = snes->type_name;
1832   PetscFunctionReturn(0);
1833 }
1834 
1835 #undef __FUNCT__
1836 #define __FUNCT__ "SNESGetSolution"
1837 /*@C
1838    SNESGetSolution - Returns the vector where the approximate solution is
1839    stored.
1840 
1841    Not Collective, but Vec is parallel if SNES is parallel
1842 
1843    Input Parameter:
1844 .  snes - the SNES context
1845 
1846    Output Parameter:
1847 .  x - the solution
1848 
1849    Level: advanced
1850 
1851 .keywords: SNES, nonlinear, get, solution
1852 
1853 .seealso: SNESGetFunction(), SNESGetSolutionUpdate()
1854 @*/
1855 int SNESGetSolution(SNES snes,Vec *x)
1856 {
1857   PetscFunctionBegin;
1858   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1859   *x = snes->vec_sol_always;
1860   PetscFunctionReturn(0);
1861 }
1862 
1863 #undef __FUNCT__
1864 #define __FUNCT__ "SNESGetSolutionUpdate"
1865 /*@C
1866    SNESGetSolutionUpdate - Returns the vector where the solution update is
1867    stored.
1868 
1869    Not Collective, but Vec is parallel if SNES is parallel
1870 
1871    Input Parameter:
1872 .  snes - the SNES context
1873 
1874    Output Parameter:
1875 .  x - the solution update
1876 
1877    Level: advanced
1878 
1879 .keywords: SNES, nonlinear, get, solution, update
1880 
1881 .seealso: SNESGetSolution(), SNESGetFunction
1882 @*/
1883 int SNESGetSolutionUpdate(SNES snes,Vec *x)
1884 {
1885   PetscFunctionBegin;
1886   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1887   *x = snes->vec_sol_update_always;
1888   PetscFunctionReturn(0);
1889 }
1890 
1891 #undef __FUNCT__
1892 #define __FUNCT__ "SNESGetFunction"
1893 /*@C
1894    SNESGetFunction - Returns the vector where the function is stored.
1895 
1896    Not Collective, but Vec is parallel if SNES is parallel
1897 
1898    Input Parameter:
1899 .  snes - the SNES context
1900 
1901    Output Parameter:
1902 +  r - the function (or PETSC_NULL)
1903 .  ctx - the function context (or PETSC_NULL)
1904 -  func - the function (or PETSC_NULL)
1905 
1906    Level: advanced
1907 
1908 .keywords: SNES, nonlinear, get, function
1909 
1910 .seealso: SNESSetFunction(), SNESGetSolution()
1911 @*/
1912 int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
1913 {
1914   PetscFunctionBegin;
1915   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1916   if (r)    *r    = snes->vec_func_always;
1917   if (ctx)  *ctx  = snes->funP;
1918   if (func) *func = snes->computefunction;
1919   PetscFunctionReturn(0);
1920 }
1921 
1922 #undef __FUNCT__
1923 #define __FUNCT__ "SNESSetOptionsPrefix"
1924 /*@C
1925    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1926    SNES options in the database.
1927 
1928    Collective on SNES
1929 
1930    Input Parameter:
1931 +  snes - the SNES context
1932 -  prefix - the prefix to prepend to all option names
1933 
1934    Notes:
1935    A hyphen (-) must NOT be given at the beginning of the prefix name.
1936    The first character of all runtime options is AUTOMATICALLY the hyphen.
1937 
1938    Level: advanced
1939 
1940 .keywords: SNES, set, options, prefix, database
1941 
1942 .seealso: SNESSetFromOptions()
1943 @*/
1944 int SNESSetOptionsPrefix(SNES snes,char *prefix)
1945 {
1946   int ierr;
1947 
1948   PetscFunctionBegin;
1949   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1950   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
1951   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1952   PetscFunctionReturn(0);
1953 }
1954 
1955 #undef __FUNCT__
1956 #define __FUNCT__ "SNESAppendOptionsPrefix"
1957 /*@C
1958    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1959    SNES options in the database.
1960 
1961    Collective on SNES
1962 
1963    Input Parameters:
1964 +  snes - the SNES context
1965 -  prefix - the prefix to prepend to all option names
1966 
1967    Notes:
1968    A hyphen (-) must NOT be given at the beginning of the prefix name.
1969    The first character of all runtime options is AUTOMATICALLY the hyphen.
1970 
1971    Level: advanced
1972 
1973 .keywords: SNES, append, options, prefix, database
1974 
1975 .seealso: SNESGetOptionsPrefix()
1976 @*/
1977 int SNESAppendOptionsPrefix(SNES snes,char *prefix)
1978 {
1979   int ierr;
1980 
1981   PetscFunctionBegin;
1982   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1983   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
1984   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
1985   PetscFunctionReturn(0);
1986 }
1987 
1988 #undef __FUNCT__
1989 #define __FUNCT__ "SNESGetOptionsPrefix"
1990 /*@C
1991    SNESGetOptionsPrefix - Sets the prefix used for searching for all
1992    SNES options in the database.
1993 
1994    Not Collective
1995 
1996    Input Parameter:
1997 .  snes - the SNES context
1998 
1999    Output Parameter:
2000 .  prefix - pointer to the prefix string used
2001 
2002    Notes: On the fortran side, the user should pass in a string 'prifix' of
2003    sufficient length to hold the prefix.
2004 
2005    Level: advanced
2006 
2007 .keywords: SNES, get, options, prefix, database
2008 
2009 .seealso: SNESAppendOptionsPrefix()
2010 @*/
2011 int SNESGetOptionsPrefix(SNES snes,char **prefix)
2012 {
2013   int ierr;
2014 
2015   PetscFunctionBegin;
2016   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2017   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2018   PetscFunctionReturn(0);
2019 }
2020 
2021 /*MC
2022    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
2023 
2024    Synopsis:
2025    int SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
2026 
2027    Not collective
2028 
2029    Input Parameters:
2030 +  name_solver - name of a new user-defined solver
2031 .  path - path (either absolute or relative) the library containing this solver
2032 .  name_create - name of routine to create method context
2033 -  routine_create - routine to create method context
2034 
2035    Notes:
2036    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
2037 
2038    If dynamic libraries are used, then the fourth input argument (routine_create)
2039    is ignored.
2040 
2041    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
2042    and others of the form ${any_environmental_variable} occuring in pathname will be
2043    replaced with appropriate values.
2044 
2045    Sample usage:
2046 .vb
2047    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2048                 "MySolverCreate",MySolverCreate);
2049 .ve
2050 
2051    Then, your solver can be chosen with the procedural interface via
2052 $     SNESSetType(snes,"my_solver")
2053    or at runtime via the option
2054 $     -snes_type my_solver
2055 
2056    Level: advanced
2057 
2058 .keywords: SNES, nonlinear, register
2059 
2060 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2061 M*/
2062 
2063 #undef __FUNCT__
2064 #define __FUNCT__ "SNESRegister"
2065 int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES))
2066 {
2067   char fullname[256];
2068   int  ierr;
2069 
2070   PetscFunctionBegin;
2071   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2072   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2073   PetscFunctionReturn(0);
2074 }
2075