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