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