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