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