xref: /petsc/src/snes/interface/snes.c (revision c38d4ed214221df9ea04de46f7761bef149d00ff)
1 /*$Id: snes.c,v 1.202 1999/11/05 14:47:05 bsmith Exp bsmith $*/
2 
3 #include "src/snes/snesimpl.h"      /*I "snes.h"  I*/
4 
5 int SNESRegisterAllCalled = 0;
6 FList SNESList = 0;
7 
8 #undef __FUNC__
9 #define __FUNC__ "SNESView"
10 /*@
11    SNESView - Prints the SNES data structure.
12 
13    Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF
14 
15    Input Parameters:
16 +  SNES - the SNES context
17 -  viewer - visualization context
18 
19    Options Database Key:
20 .  -snes_view - Calls SNESView() at end of SNESSolve()
21 
22    Notes:
23    The available visualization contexts include
24 +     VIEWER_STDOUT_SELF - standard output (default)
25 -     VIEWER_STDOUT_WORLD - synchronized standard
26          output where only the first processor opens
27          the file.  All other processors send their
28          data to the first processor to print.
29 
30    The user can open an alternative visualization context with
31    ViewerASCIIOpen() - output to a specified file.
32 
33    Level: beginner
34 
35 .keywords: SNES, view
36 
37 .seealso: ViewerASCIIOpen()
38 @*/
39 int SNESView(SNES snes,Viewer viewer)
40 {
41   SNES_KSP_EW_ConvCtx *kctx;
42   int                 ierr;
43   SLES                sles;
44   char                *type;
45   PetscTruth          isascii,isstring;
46 
47   PetscFunctionBegin;
48   PetscValidHeaderSpecific(snes,SNES_COOKIE);
49   if (!viewer) viewer = VIEWER_STDOUT_SELF;
50   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
51   PetscCheckSameComm(snes,viewer);
52 
53   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
54   ierr = PetscTypeCompare((PetscObject)viewer,STRING_VIEWER,&isstring);CHKERRQ(ierr);
55   if (isascii) {
56     ierr = ViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
57     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
58     if (type) {
59       ierr = ViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
60     } else {
61       ierr = ViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
62     }
63     if (snes->view) {
64       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
65       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
66       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
67     }
68     ierr = ViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
69     ierr = ViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
70                  snes->rtol, snes->atol, snes->trunctol, snes->xtol);CHKERRQ(ierr);
71     ierr = ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
72     ierr = ViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
73     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
74       ierr = ViewerASCIIPrintf(viewer,"  min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr);
75     }
76     if (snes->ksp_ewconv) {
77       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
78       if (kctx) {
79         ierr = ViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
80         ierr = ViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
81         ierr = ViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
82       }
83     }
84   } else if (isstring) {
85     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
86     ierr = ViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
87   }
88   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
89   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
90   ierr = SLESView(sles,viewer);CHKERRQ(ierr);
91   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 /*
96        We retain a list of functions that also take SNES command
97     line options. These are called at the end SNESSetFromOptions()
98 */
99 #define MAXSETFROMOPTIONS 5
100 static int numberofsetfromoptions;
101 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
102 
103 #undef __FUNC__
104 #define __FUNC__ "SNESAddOptionsChecker"
105 /*@
106     SNESAddOptionsChecker - Adds an additional function to check for SNES options.
107 
108     Not Collective
109 
110     Input Parameter:
111 .   snescheck - function that checks for options
112 
113     Level: developer
114 
115 .seealso: SNESSetFromOptions()
116 @*/
117 int SNESAddOptionsChecker(int (*snescheck)(SNES) )
118 {
119   PetscFunctionBegin;
120   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
121     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed");
122   }
123 
124   othersetfromoptions[numberofsetfromoptions++] = snescheck;
125   PetscFunctionReturn(0);
126 }
127 
128 #undef __FUNC__
129 #define __FUNC__ "SNESSetTypeFromOptions"
130 /*@
131    SNESSetTypeFromOptions - Sets the SNES solver type from the options database,
132         or sets a default if none is give.
133 
134    Collective on SNES
135 
136    Input Parameter:
137 .  snes - the SNES context
138 
139    Options Database Keys:
140 .  -snes_type <type> - ls, tr, umls, umtr, test
141 
142    Level: beginner
143 
144 .keywords: SNES, nonlinear, set, options, database
145 
146 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetFromOptions()
147 @*/
148 int SNESSetTypeFromOptions(SNES snes)
149 {
150   char       type[256];
151   int        ierr;
152   PetscTruth flg;
153 
154   PetscFunctionBegin;
155   PetscValidHeaderSpecific(snes,SNES_COOKIE);
156   if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp()");
157   ierr = OptionsGetString(snes->prefix,"-snes_type",type,256,&flg);CHKERRQ(ierr);
158   if (flg) {
159     ierr = SNESSetType(snes,(SNESType) type);CHKERRQ(ierr);
160   }
161   /*
162       If SNES type has not yet been set, set it now
163   */
164   if (!snes->type_name) {
165     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
166       ierr = SNESSetType(snes,SNESEQLS);CHKERRQ(ierr);
167     } else {
168       ierr = SNESSetType(snes,SNESUMTR);CHKERRQ(ierr);
169     }
170   }
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNC__
175 #define __FUNC__ "SNESSetFromOptions"
176 /*@
177    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
178 
179    Collective on SNES
180 
181    Input Parameter:
182 .  snes - the SNES context
183 
184    Options Database Keys:
185 +  -snes_type <type> - ls, tr, umls, umtr, test
186 .  -snes_stol - convergence tolerance in terms of the norm
187                 of the change in the solution between steps
188 .  -snes_atol <atol> - absolute tolerance of residual norm
189 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
190 .  -snes_max_it <max_it> - maximum number of iterations
191 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
192 .  -snes_trtol <trtol> - trust region tolerance
193 .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
194                                solver; hence iterations will continue until max_it
195                                or some other criterion is reached. Saves expense
196                                of convergence test
197 .  -snes_monitor - prints residual norm at each iteration
198 .  -snes_vecmonitor - plots solution at each iteration
199 .  -snes_vecmonitor_update - plots update to solution at each iteration
200 .  -snes_xmonitor - plots residual norm at each iteration
201 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
202 -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
203 
204     Options Database for Eisenstat-Walker method:
205 +  -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence
206 .  -snes_ksp_eq_version ver - version of  Eisenstat-Walker method
207 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
208 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
209 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
210 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
211 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
212 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
213 
214    Notes:
215    To see all options, run your program with the -help option or consult
216    the users manual.
217 
218    Level: beginner
219 
220 .keywords: SNES, nonlinear, set, options, database
221 
222 .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetTypeFromOptions()
223 @*/
224 int SNESSetFromOptions(SNES snes)
225 {
226   double     tmp;
227   SLES       sles;
228   PetscTruth flg;
229   int        ierr,i,loc[4],nmax = 4;
230   int        version   = PETSC_DEFAULT;
231   double     rtol_0    = PETSC_DEFAULT;
232   double     rtol_max  = PETSC_DEFAULT;
233   double     gamma2    = PETSC_DEFAULT;
234   double     alpha     = PETSC_DEFAULT;
235   double     alpha2    = PETSC_DEFAULT;
236   double     threshold = PETSC_DEFAULT;
237 
238   PetscFunctionBegin;
239   PetscValidHeaderSpecific(snes,SNES_COOKIE);
240   ierr = SNESSetTypeFromOptions(snes);CHKERRQ(ierr);
241 
242   loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300;
243   ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg);CHKERRQ(ierr);
244   if (flg) {
245     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
246   }
247   ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg);CHKERRQ(ierr);
248   if (flg) {
249     ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
250   }
251   ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);CHKERRQ(ierr);
252   if (flg) {
253     ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
254   }
255   ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, PETSC_NULL);CHKERRQ(ierr);
256   ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, PETSC_NULL);CHKERRQ(ierr);
257   ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg);CHKERRQ(ierr);
258   if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp);CHKERRQ(ierr); }
259   ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg);CHKERRQ(ierr);
260   if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);CHKERRQ(ierr);}
261   ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg);CHKERRQ(ierr);
262   if (flg) { snes->ksp_ewconv = 1; }
263   ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,PETSC_NULL);CHKERRQ(ierr);
264   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,PETSC_NULL);CHKERRQ(ierr);
265   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,PETSC_NULL);CHKERRQ(ierr);
266   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,PETSC_NULL);CHKERRQ(ierr);
267   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,PETSC_NULL);CHKERRQ(ierr);
268   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,PETSC_NULL);CHKERRQ(ierr);
269   ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,PETSC_NULL);CHKERRQ(ierr);
270 
271   ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg);CHKERRQ(ierr);
272   if (flg) {snes->converged = 0;}
273 
274   ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha,
275                                   alpha2,threshold);CHKERRQ(ierr);
276   ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg);CHKERRQ(ierr);
277   if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
278   ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg);CHKERRQ(ierr);
279   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
280   ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg);CHKERRQ(ierr);
281   if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
282   ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor",&flg);CHKERRQ(ierr);
283   if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
284   ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor_update",&flg);CHKERRQ(ierr);
285   if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitorUpdate,0,0);CHKERRQ(ierr);}
286   ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr);
287   if (flg) {
288     int rank;
289     ierr = MPI_Comm_rank(snes->comm,&rank);CHKERRQ(ierr);
290     if (!rank) {
291       ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
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;
1373   PetscTruth flg;
1374 
1375   PetscFunctionBegin;
1376   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1377   PetscValidHeaderSpecific(x,VEC_COOKIE);
1378   PetscCheckSameComm(snes,x);
1379   snes->vec_sol = snes->vec_sol_always = x;
1380 
1381   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);CHKERRQ(ierr);
1382   /*
1383       This version replaces the user provided Jacobian matrix with a
1384       matrix-free version but still employs the user-provided preconditioner matrix
1385   */
1386   if (flg) {
1387     Mat J;
1388     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1389     PLogObjectParent(snes,J);
1390     snes->mfshell = J;
1391     snes->jacobian = J;
1392     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1393       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1394     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1395       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1396     } else {
1397       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option");
1398     }
1399     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1400   }
1401   ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);CHKERRQ(ierr);
1402   /*
1403       This version replaces both the user-provided Jacobian and the user-
1404       provided preconditioner matrix with the default matrix free version.
1405    */
1406   if (flg) {
1407     Mat J;
1408     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1409     PLogObjectParent(snes,J);
1410     snes->mfshell = J;
1411     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1412       ierr = SNESSetJacobian(snes,J,J,0,snes->funP);CHKERRQ(ierr);
1413       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1414     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1415       ierr = SNESSetHessian(snes,J,J,0,snes->funP);CHKERRQ(ierr);
1416       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1417     } else {
1418       SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option");
1419     }
1420     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1421   }
1422   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
1423     PetscTruth iseqtr;
1424 
1425     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1426     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first");
1427     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1428     if (snes->vec_func == snes->vec_sol) {
1429       SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector");
1430     }
1431 
1432     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
1433     ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr);
1434     if (snes->ksp_ewconv && !iseqtr) {
1435       SLES sles; KSP ksp;
1436       ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1437       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1438       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,(void *)snes);CHKERRQ(ierr);
1439     }
1440   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1441     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1442     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first");
1443     if (!snes->computeumfunction) {
1444       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first");
1445     }
1446     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian()");
1447   } else {
1448     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class");
1449   }
1450   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
1451   snes->setupcalled = 1;
1452   PetscFunctionReturn(0);
1453 }
1454 
1455 #undef __FUNC__
1456 #define __FUNC__ "SNESDestroy"
1457 /*@C
1458    SNESDestroy - Destroys the nonlinear solver context that was created
1459    with SNESCreate().
1460 
1461    Collective on SNES
1462 
1463    Input Parameter:
1464 .  snes - the SNES context
1465 
1466    Level: beginner
1467 
1468 .keywords: SNES, nonlinear, destroy
1469 
1470 .seealso: SNESCreate(), SNESSolve()
1471 @*/
1472 int SNESDestroy(SNES snes)
1473 {
1474   int i,ierr;
1475 
1476   PetscFunctionBegin;
1477   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1478   if (--snes->refct > 0) PetscFunctionReturn(0);
1479 
1480   /* if memory was published with AMS then destroy it */
1481   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1482 
1483   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1484   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
1485   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
1486   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1487   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1488   for (i=0; i<snes->numbermonitors; i++ ) {
1489     if (snes->monitordestroy[i]) {
1490       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1491     }
1492   }
1493   PLogObjectDestroy((PetscObject)snes);
1494   PetscHeaderDestroy((PetscObject)snes);
1495   PetscFunctionReturn(0);
1496 }
1497 
1498 /* ----------- Routines to set solver parameters ---------- */
1499 
1500 #undef __FUNC__
1501 #define __FUNC__ "SNESSetTolerances"
1502 /*@
1503    SNESSetTolerances - Sets various parameters used in convergence tests.
1504 
1505    Collective on SNES
1506 
1507    Input Parameters:
1508 +  snes - the SNES context
1509 .  atol - absolute convergence tolerance
1510 .  rtol - relative convergence tolerance
1511 .  stol -  convergence tolerance in terms of the norm
1512            of the change in the solution between steps
1513 .  maxit - maximum number of iterations
1514 -  maxf - maximum number of function evaluations
1515 
1516    Options Database Keys:
1517 +    -snes_atol <atol> - Sets atol
1518 .    -snes_rtol <rtol> - Sets rtol
1519 .    -snes_stol <stol> - Sets stol
1520 .    -snes_max_it <maxit> - Sets maxit
1521 -    -snes_max_funcs <maxf> - Sets maxf
1522 
1523    Notes:
1524    The default maximum number of iterations is 50.
1525    The default maximum number of function evaluations is 1000.
1526 
1527    Level: intermediate
1528 
1529 .keywords: SNES, nonlinear, set, convergence, tolerances
1530 
1531 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1532 @*/
1533 int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf)
1534 {
1535   PetscFunctionBegin;
1536   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1537   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1538   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1539   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1540   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1541   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1542   PetscFunctionReturn(0);
1543 }
1544 
1545 #undef __FUNC__
1546 #define __FUNC__ "SNESGetTolerances"
1547 /*@
1548    SNESGetTolerances - Gets various parameters used in convergence tests.
1549 
1550    Not Collective
1551 
1552    Input Parameters:
1553 +  snes - the SNES context
1554 .  atol - absolute convergence tolerance
1555 .  rtol - relative convergence tolerance
1556 .  stol -  convergence tolerance in terms of the norm
1557            of the change in the solution between steps
1558 .  maxit - maximum number of iterations
1559 -  maxf - maximum number of function evaluations
1560 
1561    Notes:
1562    The user can specify PETSC_NULL for any parameter that is not needed.
1563 
1564    Level: intermediate
1565 
1566 .keywords: SNES, nonlinear, get, convergence, tolerances
1567 
1568 .seealso: SNESSetTolerances()
1569 @*/
1570 int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf)
1571 {
1572   PetscFunctionBegin;
1573   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1574   if (atol)  *atol  = snes->atol;
1575   if (rtol)  *rtol  = snes->rtol;
1576   if (stol)  *stol  = snes->xtol;
1577   if (maxit) *maxit = snes->max_its;
1578   if (maxf)  *maxf  = snes->max_funcs;
1579   PetscFunctionReturn(0);
1580 }
1581 
1582 #undef __FUNC__
1583 #define __FUNC__ "SNESSetTrustRegionTolerance"
1584 /*@
1585    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1586 
1587    Collective on SNES
1588 
1589    Input Parameters:
1590 +  snes - the SNES context
1591 -  tol - tolerance
1592 
1593    Options Database Key:
1594 .  -snes_trtol <tol> - Sets tol
1595 
1596    Level: intermediate
1597 
1598 .keywords: SNES, nonlinear, set, trust region, tolerance
1599 
1600 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1601 @*/
1602 int SNESSetTrustRegionTolerance(SNES snes,double tol)
1603 {
1604   PetscFunctionBegin;
1605   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1606   snes->deltatol = tol;
1607   PetscFunctionReturn(0);
1608 }
1609 
1610 #undef __FUNC__
1611 #define __FUNC__ "SNESSetMinimizationFunctionTolerance"
1612 /*@
1613    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
1614    for unconstrained minimization solvers.
1615 
1616    Collective on SNES
1617 
1618    Input Parameters:
1619 +  snes - the SNES context
1620 -  ftol - minimum function tolerance
1621 
1622    Options Database Key:
1623 .  -snes_fmin <ftol> - Sets ftol
1624 
1625    Note:
1626    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
1627    methods only.
1628 
1629    Level: intermediate
1630 
1631 .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
1632 
1633 .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
1634 @*/
1635 int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol)
1636 {
1637   PetscFunctionBegin;
1638   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1639   snes->fmin = ftol;
1640   PetscFunctionReturn(0);
1641 }
1642 /*
1643    Duplicate the lg monitors for SNES from KSP; for some reason with
1644    dynamic libraries things don't work under Sun4 if we just use
1645    macros instead of functions
1646 */
1647 #undef __FUNC__
1648 #define __FUNC__ "SNESLGMonitor"
1649 int SNESLGMonitor(SNES snes,int it,double norm,void *ctx)
1650 {
1651   int ierr;
1652 
1653   PetscFunctionBegin;
1654   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1655   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1656   PetscFunctionReturn(0);
1657 }
1658 
1659 #undef __FUNC__
1660 #define __FUNC__ "SNESLGMonitorCreate"
1661 int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n, DrawLG *draw)
1662 {
1663   int ierr;
1664 
1665   PetscFunctionBegin;
1666   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1667   PetscFunctionReturn(0);
1668 }
1669 
1670 #undef __FUNC__
1671 #define __FUNC__ "SNESLGMonitorDestroy"
1672 int SNESLGMonitorDestroy(DrawLG draw)
1673 {
1674   int ierr;
1675 
1676   PetscFunctionBegin;
1677   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1678   PetscFunctionReturn(0);
1679 }
1680 
1681 /* ------------ Routines to set performance monitoring options ----------- */
1682 
1683 #undef __FUNC__
1684 #define __FUNC__ "SNESSetMonitor"
1685 /*@C
1686    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
1687    iteration of the nonlinear solver to display the iteration's
1688    progress.
1689 
1690    Collective on SNES
1691 
1692    Input Parameters:
1693 +  snes - the SNES context
1694 .  func - monitoring routine
1695 .  mctx - [optional] user-defined context for private data for the
1696           monitor routine (may be PETSC_NULL)
1697 -  monitordestroy - options routine that frees monitor context
1698 
1699    Calling sequence of func:
1700 $     int func(SNES snes,int its, double norm,void *mctx)
1701 
1702 +    snes - the SNES context
1703 .    its - iteration number
1704 .    norm - 2-norm function value (may be estimated)
1705             for SNES_NONLINEAR_EQUATIONS methods
1706 .    norm - 2-norm gradient value (may be estimated)
1707             for SNES_UNCONSTRAINED_MINIMIZATION methods
1708 -    mctx - [optional] monitoring context
1709 
1710    Options Database Keys:
1711 +    -snes_monitor        - sets SNESDefaultMonitor()
1712 .    -snes_xmonitor       - sets line graph monitor,
1713                             uses SNESLGMonitorCreate()
1714 _    -snes_cancelmonitors - cancels all monitors that have
1715                             been hardwired into a code by
1716                             calls to SNESSetMonitor(), but
1717                             does not cancel those set via
1718                             the options database.
1719 
1720    Notes:
1721    Several different monitoring routines may be set by calling
1722    SNESSetMonitor() multiple times; all will be called in the
1723    order in which they were set.
1724 
1725    Level: intermediate
1726 
1727 .keywords: SNES, nonlinear, set, monitor
1728 
1729 .seealso: SNESDefaultMonitor(), SNESClearMonitor()
1730 @*/
1731 int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx,int (*monitordestroy)(void *))
1732 {
1733   PetscFunctionBegin;
1734   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1735   if (snes->numbermonitors >= MAXSNESMONITORS) {
1736     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set");
1737   }
1738 
1739   snes->monitor[snes->numbermonitors]           = func;
1740   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1741   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1742   PetscFunctionReturn(0);
1743 }
1744 
1745 #undef __FUNC__
1746 #define __FUNC__ "SNESClearMonitor"
1747 /*@C
1748    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1749 
1750    Collective on SNES
1751 
1752    Input Parameters:
1753 .  snes - the SNES context
1754 
1755    Options Database:
1756 .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1757     into a code by calls to SNESSetMonitor(), but does not cancel those
1758     set via the options database
1759 
1760    Notes:
1761    There is no way to clear one specific monitor from a SNES object.
1762 
1763    Level: intermediate
1764 
1765 .keywords: SNES, nonlinear, set, monitor
1766 
1767 .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1768 @*/
1769 int SNESClearMonitor( SNES snes )
1770 {
1771   PetscFunctionBegin;
1772   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1773   snes->numbermonitors = 0;
1774   PetscFunctionReturn(0);
1775 }
1776 
1777 #undef __FUNC__
1778 #define __FUNC__ "SNESSetConvergenceTest"
1779 /*@C
1780    SNESSetConvergenceTest - Sets the function that is to be used
1781    to test for convergence of the nonlinear iterative solution.
1782 
1783    Collective on SNES
1784 
1785    Input Parameters:
1786 +  snes - the SNES context
1787 .  func - routine to test for convergence
1788 -  cctx - [optional] context for private data for the convergence routine
1789           (may be PETSC_NULL)
1790 
1791    Calling sequence of func:
1792 $     int func (SNES snes,double xnorm,double gnorm,double f,SNESConvergedReason *reason,void *cctx)
1793 
1794 +    snes - the SNES context
1795 .    cctx - [optional] convergence context
1796 .    reason - reason for convergence/divergence
1797 .    xnorm - 2-norm of current iterate
1798 .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1799 .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1800 .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1801 -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
1802 
1803    Level: advanced
1804 
1805 .keywords: SNES, nonlinear, set, convergence, test
1806 
1807 .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1808           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1809 @*/
1810 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,SNESConvergedReason*,void*),void *cctx)
1811 {
1812   PetscFunctionBegin;
1813   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1814   (snes)->converged = func;
1815   (snes)->cnvP      = cctx;
1816   PetscFunctionReturn(0);
1817 }
1818 
1819 #undef __FUNC__
1820 #define __FUNC__ "SNESGetConvergedReason"
1821 /*@C
1822    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1823 
1824    Not Collective
1825 
1826    Input Parameter:
1827 .  snes - the SNES context
1828 
1829    Output Parameter:
1830 .  reason - negative value indicates diverged, positive value converged, see snes.h or the
1831             manual pages for the individual convergence tests for complete lists
1832 
1833    Level: intermediate
1834 
1835    Notes: Can only be called after the call the SNESSolve() is complete.
1836 
1837 .keywords: SNES, nonlinear, set, convergence, test
1838 
1839 .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1840           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1841 @*/
1842 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1843 {
1844   PetscFunctionBegin;
1845   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1846   *reason = snes->reason;
1847   PetscFunctionReturn(0);
1848 }
1849 
1850 #undef __FUNC__
1851 #define __FUNC__ "SNESSetConvergenceHistory"
1852 /*@
1853    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1854 
1855    Collective on SNES
1856 
1857    Input Parameters:
1858 +  snes - iterative context obtained from SNESCreate()
1859 .  a   - array to hold history
1860 .  its - integer array holds the number of linear iterations (or
1861          negative if not converged) for each solve.
1862 .  na  - size of a and its
1863 -  reset - PETSC_TRUTH indicates each new nonlinear solve resets the history counter to zero,
1864            else it continues storing new values for new nonlinear solves after the old ones
1865 
1866    Notes:
1867    If set, this array will contain the function norms (for
1868    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1869    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1870    at each step.
1871 
1872    This routine is useful, e.g., when running a code for purposes
1873    of accurate performance monitoring, when no I/O should be done
1874    during the section of code that is being timed.
1875 
1876    Level: intermediate
1877 
1878 .keywords: SNES, set, convergence, history
1879 
1880 .seealso: SNESGetConvergenceHistory()
1881 
1882 @*/
1883 int SNESSetConvergenceHistory(SNES snes, double *a, int *its,int na,PetscTruth reset)
1884 {
1885   PetscFunctionBegin;
1886   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1887   if (na) PetscValidScalarPointer(a);
1888   snes->conv_hist       = a;
1889   snes->conv_hist_its   = its;
1890   snes->conv_hist_max   = na;
1891   snes->conv_hist_reset = reset;
1892   PetscFunctionReturn(0);
1893 }
1894 
1895 #undef __FUNC__
1896 #define __FUNC__ "SNESGetConvergenceHistory"
1897 /*@C
1898    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1899 
1900    Collective on SNES
1901 
1902    Input Parameter:
1903 .  snes - iterative context obtained from SNESCreate()
1904 
1905    Output Parameters:
1906 .  a   - array to hold history
1907 .  its - integer array holds the number of linear iterations (or
1908          negative if not converged) for each solve.
1909 -  na  - size of a and its
1910 
1911    Notes:
1912     The calling sequence for this routine in Fortran is
1913 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1914 
1915    This routine is useful, e.g., when running a code for purposes
1916    of accurate performance monitoring, when no I/O should be done
1917    during the section of code that is being timed.
1918 
1919    Level: intermediate
1920 
1921 .keywords: SNES, get, convergence, history
1922 
1923 .seealso: SNESSetConvergencHistory()
1924 
1925 @*/
1926 int SNESGetConvergenceHistory(SNES snes, double **a, int **its,int *na)
1927 {
1928   PetscFunctionBegin;
1929   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1930   if (a)   *a   = snes->conv_hist;
1931   if (its) *its = snes->conv_hist_its;
1932   if (na) *na   = snes->conv_hist_len;
1933   PetscFunctionReturn(0);
1934 }
1935 
1936 #undef __FUNC__
1937 #define __FUNC__ "SNESScaleStep_Private"
1938 /*
1939    SNESScaleStep_Private - Scales a step so that its length is less than the
1940    positive parameter delta.
1941 
1942     Input Parameters:
1943 +   snes - the SNES context
1944 .   y - approximate solution of linear system
1945 .   fnorm - 2-norm of current function
1946 -   delta - trust region size
1947 
1948     Output Parameters:
1949 +   gpnorm - predicted function norm at the new point, assuming local
1950     linearization.  The value is zero if the step lies within the trust
1951     region, and exceeds zero otherwise.
1952 -   ynorm - 2-norm of the step
1953 
1954     Note:
1955     For non-trust region methods such as SNESEQLS, the parameter delta
1956     is set to be the maximum allowable step size.
1957 
1958 .keywords: SNES, nonlinear, scale, step
1959 */
1960 int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta,
1961                           double *gpnorm,double *ynorm)
1962 {
1963   double norm;
1964   Scalar cnorm;
1965   int    ierr;
1966 
1967   PetscFunctionBegin;
1968   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1969   PetscValidHeaderSpecific(y,VEC_COOKIE);
1970   PetscCheckSameComm(snes,y);
1971 
1972   ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr);
1973   if (norm > *delta) {
1974      norm = *delta/norm;
1975      *gpnorm = (1.0 - norm)*(*fnorm);
1976      cnorm = norm;
1977      VecScale( &cnorm, y );
1978      *ynorm = *delta;
1979   } else {
1980      *gpnorm = 0.0;
1981      *ynorm = norm;
1982   }
1983   PetscFunctionReturn(0);
1984 }
1985 
1986 #undef __FUNC__
1987 #define __FUNC__ "SNESSolve"
1988 /*@
1989    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1990    SNESCreate() and optional routines of the form SNESSetXXX().
1991 
1992    Collective on SNES
1993 
1994    Input Parameters:
1995 +  snes - the SNES context
1996 -  x - the solution vector
1997 
1998    Output Parameter:
1999 .  its - number of iterations until termination
2000 
2001    Notes:
2002    The user should initialize the vector, x, with the initial guess
2003    for the nonlinear solve prior to calling SNESSolve.  In particular,
2004    to employ an initial guess of zero, the user should explicitly set
2005    this vector to zero by calling VecSet().
2006 
2007    Level: beginner
2008 
2009 .keywords: SNES, nonlinear, solve
2010 
2011 .seealso: SNESCreate(), SNESDestroy()
2012 @*/
2013 int SNESSolve(SNES snes,Vec x,int *its)
2014 {
2015   int        ierr;
2016   PetscTruth flg;
2017 
2018   PetscFunctionBegin;
2019   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2020   PetscValidHeaderSpecific(x,VEC_COOKIE);
2021   PetscCheckSameComm(snes,x);
2022   PetscValidIntPointer(its);
2023   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
2024   else {snes->vec_sol = snes->vec_sol_always = x;}
2025   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
2026   PLogEventBegin(SNES_Solve,snes,0,0,0);
2027   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
2028   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
2029   PLogEventEnd(SNES_Solve,snes,0,0,0);
2030   ierr = OptionsHasName(snes->prefix,"-snes_view", &flg);CHKERRQ(ierr);
2031   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD);CHKERRQ(ierr); }
2032   PetscFunctionReturn(0);
2033 }
2034 
2035 /* --------- Internal routines for SNES Package --------- */
2036 
2037 #undef __FUNC__
2038 #define __FUNC__ "SNESSetType"
2039 /*@C
2040    SNESSetType - Sets the method for the nonlinear solver.
2041 
2042    Collective on SNES
2043 
2044    Input Parameters:
2045 +  snes - the SNES context
2046 -  type - a known method
2047 
2048    Options Database Key:
2049 .  -snes_type <type> - Sets the method; use -help for a list
2050    of available methods (for instance, ls or tr)
2051 
2052    Notes:
2053    See "petsc/include/snes.h" for available methods (for instance)
2054 +    SNESEQLS - Newton's method with line search
2055      (systems of nonlinear equations)
2056 .    SNESEQTR - Newton's method with trust region
2057      (systems of nonlinear equations)
2058 .    SNESUMTR - Newton's method with trust region
2059      (unconstrained minimization)
2060 -    SNESUMLS - Newton's method with line search
2061      (unconstrained minimization)
2062 
2063   Normally, it is best to use the SNESSetFromOptions() command and then
2064   set the SNES solver type from the options database rather than by using
2065   this routine.  Using the options database provides the user with
2066   maximum flexibility in evaluating the many nonlinear solvers.
2067   The SNESSetType() routine is provided for those situations where it
2068   is necessary to set the nonlinear solver independently of the command
2069   line or options database.  This might be the case, for example, when
2070   the choice of solver changes during the execution of the program,
2071   and the user's application is taking responsibility for choosing the
2072   appropriate method.  In other words, this routine is not for beginners.
2073 
2074   Level: intermediate
2075 
2076 .keywords: SNES, set, type
2077 @*/
2078 int SNESSetType(SNES snes,SNESType type)
2079 {
2080   int        ierr, (*r)(SNES);
2081   PetscTruth match;
2082 
2083   PetscFunctionBegin;
2084   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2085   PetscValidCharPointer(type);
2086 
2087   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2088   if (match) PetscFunctionReturn(0);
2089 
2090   if (snes->setupcalled) {
2091     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
2092     snes->data = 0;
2093   }
2094 
2095   /* Get the function pointers for the iterative method requested */
2096   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2097 
2098   ierr =  FListFind(snes->comm, SNESList, type,(int (**)(void *)) &r );CHKERRQ(ierr);
2099 
2100   if (!r) SETERRQ1(1,1,"Unable to find requested SNES type %s",type);
2101 
2102   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
2103   snes->data = 0;
2104   ierr = (*r)(snes);CHKERRQ(ierr);
2105 
2106   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2107   snes->set_method_called = 1;
2108 
2109   PetscFunctionReturn(0);
2110 }
2111 
2112 
2113 /* --------------------------------------------------------------------- */
2114 #undef __FUNC__
2115 #define __FUNC__ "SNESRegisterDestroy"
2116 /*@C
2117    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2118    registered by SNESRegisterDynamic().
2119 
2120    Not Collective
2121 
2122    Level: advanced
2123 
2124 .keywords: SNES, nonlinear, register, destroy
2125 
2126 .seealso: SNESRegisterAll(), SNESRegisterAll()
2127 @*/
2128 int SNESRegisterDestroy(void)
2129 {
2130   int ierr;
2131 
2132   PetscFunctionBegin;
2133   if (SNESList) {
2134     ierr = FListDestroy( SNESList );CHKERRQ(ierr);
2135     SNESList = 0;
2136   }
2137   SNESRegisterAllCalled = 0;
2138   PetscFunctionReturn(0);
2139 }
2140 
2141 #undef __FUNC__
2142 #define __FUNC__ "SNESGetType"
2143 /*@C
2144    SNESGetType - Gets the SNES method type and name (as a string).
2145 
2146    Not Collective
2147 
2148    Input Parameter:
2149 .  snes - nonlinear solver context
2150 
2151    Output Parameter:
2152 .  type - SNES method (a charactor string)
2153 
2154    Level: intermediate
2155 
2156 .keywords: SNES, nonlinear, get, type, name
2157 @*/
2158 int SNESGetType(SNES snes, SNESType *type)
2159 {
2160   PetscFunctionBegin;
2161   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2162   *type = snes->type_name;
2163   PetscFunctionReturn(0);
2164 }
2165 
2166 #undef __FUNC__
2167 #define __FUNC__ "SNESGetSolution"
2168 /*@C
2169    SNESGetSolution - Returns the vector where the approximate solution is
2170    stored.
2171 
2172    Not Collective, but Vec is parallel if SNES is parallel
2173 
2174    Input Parameter:
2175 .  snes - the SNES context
2176 
2177    Output Parameter:
2178 .  x - the solution
2179 
2180    Level: advanced
2181 
2182 .keywords: SNES, nonlinear, get, solution
2183 
2184 .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
2185 @*/
2186 int SNESGetSolution(SNES snes,Vec *x)
2187 {
2188   PetscFunctionBegin;
2189   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2190   *x = snes->vec_sol_always;
2191   PetscFunctionReturn(0);
2192 }
2193 
2194 #undef __FUNC__
2195 #define __FUNC__ "SNESGetSolutionUpdate"
2196 /*@C
2197    SNESGetSolutionUpdate - Returns the vector where the solution update is
2198    stored.
2199 
2200    Not Collective, but Vec is parallel if SNES is parallel
2201 
2202    Input Parameter:
2203 .  snes - the SNES context
2204 
2205    Output Parameter:
2206 .  x - the solution update
2207 
2208    Level: advanced
2209 
2210 .keywords: SNES, nonlinear, get, solution, update
2211 
2212 .seealso: SNESGetSolution(), SNESGetFunction
2213 @*/
2214 int SNESGetSolutionUpdate(SNES snes,Vec *x)
2215 {
2216   PetscFunctionBegin;
2217   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2218   *x = snes->vec_sol_update_always;
2219   PetscFunctionReturn(0);
2220 }
2221 
2222 #undef __FUNC__
2223 #define __FUNC__ "SNESGetFunction"
2224 /*@C
2225    SNESGetFunction - Returns the vector where the function is stored.
2226 
2227    Not Collective, but Vec is parallel if SNES is parallel
2228 
2229    Input Parameter:
2230 .  snes - the SNES context
2231 
2232    Output Parameter:
2233 +  r - the function (or PETSC_NULL)
2234 -  ctx - the function context (or PETSC_NULL)
2235 
2236    Notes:
2237    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
2238    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
2239    SNESGetMinimizationFunction() and SNESGetGradient();
2240 
2241    Level: advanced
2242 
2243 .keywords: SNES, nonlinear, get, function
2244 
2245 .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
2246           SNESGetGradient()
2247 
2248 @*/
2249 int SNESGetFunction(SNES snes,Vec *r,void **ctx)
2250 {
2251   PetscFunctionBegin;
2252   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2253   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
2254     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
2255   }
2256   if (r)   *r = snes->vec_func_always;
2257   if (ctx) *ctx = snes->funP;
2258   PetscFunctionReturn(0);
2259 }
2260 
2261 #undef __FUNC__
2262 #define __FUNC__ "SNESGetGradient"
2263 /*@C
2264    SNESGetGradient - Returns the vector where the gradient is stored.
2265 
2266    Not Collective, but Vec is parallel if SNES is parallel
2267 
2268    Input Parameter:
2269 .  snes - the SNES context
2270 
2271    Output Parameter:
2272 +  r - the gradient (or PETSC_NULL)
2273 -  ctx - the gradient context (or PETSC_NULL)
2274 
2275    Notes:
2276    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
2277    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
2278    SNESGetFunction().
2279 
2280    Level: advanced
2281 
2282 .keywords: SNES, nonlinear, get, gradient
2283 
2284 .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(),
2285           SNESSetGradient(), SNESSetFunction()
2286 
2287 @*/
2288 int SNESGetGradient(SNES snes,Vec *r,void **ctx)
2289 {
2290   PetscFunctionBegin;
2291   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2292   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2293     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
2294   }
2295   if (r)   *r = snes->vec_func_always;
2296   if (ctx) *ctx = snes->funP;
2297   PetscFunctionReturn(0);
2298 }
2299 
2300 #undef __FUNC__
2301 #define __FUNC__ "SNESGetMinimizationFunction"
2302 /*@C
2303    SNESGetMinimizationFunction - Returns the scalar function value for
2304    unconstrained minimization problems.
2305 
2306    Not Collective
2307 
2308    Input Parameter:
2309 .  snes - the SNES context
2310 
2311    Output Parameter:
2312 +  r - the function (or PETSC_NULL)
2313 -  ctx - the function context (or PETSC_NULL)
2314 
2315    Notes:
2316    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
2317    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
2318    SNESGetFunction().
2319 
2320    Level: advanced
2321 
2322 .keywords: SNES, nonlinear, get, function
2323 
2324 .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction()
2325 
2326 @*/
2327 int SNESGetMinimizationFunction(SNES snes,double *r,void **ctx)
2328 {
2329   PetscFunctionBegin;
2330   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2331   PetscValidScalarPointer(r);
2332   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2333     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only");
2334   }
2335   if (r)   *r = snes->fc;
2336   if (ctx) *ctx = snes->umfunP;
2337   PetscFunctionReturn(0);
2338 }
2339 
2340 #undef __FUNC__
2341 #define __FUNC__ "SNESSetOptionsPrefix"
2342 /*@C
2343    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2344    SNES options in the database.
2345 
2346    Collective on SNES
2347 
2348    Input Parameter:
2349 +  snes - the SNES context
2350 -  prefix - the prefix to prepend to all option names
2351 
2352    Notes:
2353    A hyphen (-) must NOT be given at the beginning of the prefix name.
2354    The first character of all runtime options is AUTOMATICALLY the hyphen.
2355 
2356    Level: advanced
2357 
2358 .keywords: SNES, set, options, prefix, database
2359 
2360 .seealso: SNESSetFromOptions()
2361 @*/
2362 int SNESSetOptionsPrefix(SNES snes,char *prefix)
2363 {
2364   int ierr;
2365 
2366   PetscFunctionBegin;
2367   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2368   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr);
2369   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2370   PetscFunctionReturn(0);
2371 }
2372 
2373 #undef __FUNC__
2374 #define __FUNC__ "SNESAppendOptionsPrefix"
2375 /*@C
2376    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2377    SNES options in the database.
2378 
2379    Collective on SNES
2380 
2381    Input Parameters:
2382 +  snes - the SNES context
2383 -  prefix - the prefix to prepend to all option names
2384 
2385    Notes:
2386    A hyphen (-) must NOT be given at the beginning of the prefix name.
2387    The first character of all runtime options is AUTOMATICALLY the hyphen.
2388 
2389    Level: advanced
2390 
2391 .keywords: SNES, append, options, prefix, database
2392 
2393 .seealso: SNESGetOptionsPrefix()
2394 @*/
2395 int SNESAppendOptionsPrefix(SNES snes,char *prefix)
2396 {
2397   int ierr;
2398 
2399   PetscFunctionBegin;
2400   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2401   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr);
2402   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
2403   PetscFunctionReturn(0);
2404 }
2405 
2406 #undef __FUNC__
2407 #define __FUNC__ "SNESGetOptionsPrefix"
2408 /*@C
2409    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2410    SNES options in the database.
2411 
2412    Not Collective
2413 
2414    Input Parameter:
2415 .  snes - the SNES context
2416 
2417    Output Parameter:
2418 .  prefix - pointer to the prefix string used
2419 
2420    Notes: On the fortran side, the user should pass in a string 'prifix' of
2421    sufficient length to hold the prefix.
2422 
2423    Level: advanced
2424 
2425 .keywords: SNES, get, options, prefix, database
2426 
2427 .seealso: SNESAppendOptionsPrefix()
2428 @*/
2429 int SNESGetOptionsPrefix(SNES snes,char **prefix)
2430 {
2431   int ierr;
2432 
2433   PetscFunctionBegin;
2434   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2435   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr);
2436   PetscFunctionReturn(0);
2437 }
2438 
2439 #undef __FUNC__
2440 #define __FUNC__ "SNESPrintHelp"
2441 /*@
2442    SNESPrintHelp - Prints all options for the SNES component.
2443 
2444    Collective on SNES
2445 
2446    Input Parameter:
2447 .  snes - the SNES context
2448 
2449    Options Database Keys:
2450 +  -help - Prints SNES options
2451 -  -h - Prints SNES options
2452 
2453    Level: beginner
2454 
2455 .keywords: SNES, nonlinear, help
2456 
2457 .seealso: SNESSetFromOptions()
2458 @*/
2459 int SNESPrintHelp(SNES snes)
2460 {
2461   char                p[64];
2462   SNES_KSP_EW_ConvCtx *kctx;
2463   int                 ierr;
2464 
2465   PetscFunctionBegin;
2466   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2467 
2468   ierr = PetscStrcpy(p,"-");CHKERRQ(ierr);
2469   if (snes->prefix) PetscStrcat(p, snes->prefix);
2470 
2471   kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
2472 
2473   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2474   ierr = (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");CHKERRQ(ierr);
2475   ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr);
2476   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);CHKERRQ(ierr);
2477   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);CHKERRQ(ierr);
2478   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);CHKERRQ(ierr);
2479   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);CHKERRQ(ierr);
2480   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);CHKERRQ(ierr);
2481   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);CHKERRQ(ierr);
2482   ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);CHKERRQ(ierr);
2483   ierr = (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");CHKERRQ(ierr);
2484   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);CHKERRQ(ierr);
2485   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_monitor: use default SNES convergence monitor, prints\n\
2486     residual norm at each iteration.\n",p);CHKERRQ(ierr);
2487   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_smonitor: same as the above, but prints fewer digits of the\n\
2488     residual norm for small residual norms. This is useful to conceal\n\
2489     meaningless digits that may be different on different machines.\n",p);CHKERRQ(ierr);
2490   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);CHKERRQ(ierr);
2491   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_vecmonitor: plots solution at each iteration \n",p);CHKERRQ(ierr);
2492   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_vecmonitor_update: plots update to solution at each iteration \n",p);CHKERRQ(ierr);
2493   if (snes->type == SNES_NONLINEAR_EQUATIONS) {
2494     ierr = (*PetscHelpPrintf)(snes->comm,
2495      " Options for solving systems of nonlinear equations only:\n");CHKERRQ(ierr);
2496     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Jacobian\n",p);CHKERRQ(ierr);
2497     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Jacobian\n",p);CHKERRQ(ierr);
2498     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);CHKERRQ(ierr);
2499     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_ksp_monitor - if using matrix-free multiply then prints h at each KSP iteration\n",p);CHKERRQ(ierr);
2500     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);CHKERRQ(ierr);
2501     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);CHKERRQ(ierr);
2502     ierr = (*PetscHelpPrintf)(snes->comm,
2503      "     %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);CHKERRQ(ierr);
2504     ierr = (*PetscHelpPrintf)(snes->comm,
2505      "     %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);CHKERRQ(ierr);
2506     ierr = (*PetscHelpPrintf)(snes->comm,
2507      "     %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);CHKERRQ(ierr);
2508     ierr = (*PetscHelpPrintf)(snes->comm,
2509      "     %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);CHKERRQ(ierr);
2510     ierr = (*PetscHelpPrintf)(snes->comm,
2511      "     %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);CHKERRQ(ierr);
2512     ierr = (*PetscHelpPrintf)(snes->comm,
2513      "     %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);CHKERRQ(ierr);
2514     ierr = (*PetscHelpPrintf)(snes->comm,
2515      "     %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);CHKERRQ(ierr);
2516   } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) {
2517     ierr = (*PetscHelpPrintf)(snes->comm," Options for solving unconstrained minimization problems only:\n");CHKERRQ(ierr);
2518     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);CHKERRQ(ierr);
2519     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_fd: use finite differences for Hessian\n",p);CHKERRQ(ierr);
2520     ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf: use matrix-free Hessian\n",p);CHKERRQ(ierr);
2521   }
2522   ierr = (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <type> for help on ",p);CHKERRQ(ierr);
2523   ierr = (*PetscHelpPrintf)(snes->comm,"a particular method\n");CHKERRQ(ierr);
2524   if (snes->printhelp) {
2525     ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr);
2526   }
2527   PetscFunctionReturn(0);
2528 }
2529 
2530 /*MC
2531    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
2532 
2533    Synopsis:
2534    SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
2535 
2536    Not collective
2537 
2538    Input Parameters:
2539 +  name_solver - name of a new user-defined solver
2540 .  path - path (either absolute or relative) the library containing this solver
2541 .  name_create - name of routine to create method context
2542 -  routine_create - routine to create method context
2543 
2544    Notes:
2545    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
2546 
2547    If dynamic libraries are used, then the fourth input argument (routine_create)
2548    is ignored.
2549 
2550    Sample usage:
2551 .vb
2552    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2553                 "MySolverCreate",MySolverCreate);
2554 .ve
2555 
2556    Then, your solver can be chosen with the procedural interface via
2557 $     SNESSetType(snes,"my_solver")
2558    or at runtime via the option
2559 $     -snes_type my_solver
2560 
2561    Level: advanced
2562 
2563    $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values.
2564 
2565 .keywords: SNES, nonlinear, register
2566 
2567 .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2568 M*/
2569 
2570 #undef __FUNC__
2571 #define __FUNC__ "SNESRegister"
2572 int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES))
2573 {
2574   char fullname[256];
2575   int  ierr;
2576 
2577   PetscFunctionBegin;
2578   ierr = FListConcat(path,name,fullname); CHKERRQ(ierr);
2579   ierr = FListAdd(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr);
2580   PetscFunctionReturn(0);
2581 }
2582