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