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