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