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