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