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