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