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