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