xref: /petsc/src/snes/interface/snes.c (revision 4482741e5b2e2bbc854fb1f8dba65221386520f2)
1 /*$Id: snes.c,v 1.235 2001/08/21 21:03:49 bsmith Exp $*/
2 
3 #include "src/snes/snesimpl.h"      /*I "petscsnes.h"  I*/
4 
5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6 PetscFList SNESList              = PETSC_NULL;
7 
8 /* Logging support */
9 int SNES_COOKIE = 0;
10 int SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0;
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "SNESView"
14 /*@C
15    SNESView - Prints the SNES data structure.
16 
17    Collective on SNES
18 
19    Input Parameters:
20 +  SNES - the SNES context
21 -  viewer - visualization context
22 
23    Options Database Key:
24 .  -snes_view - Calls SNESView() at end of SNESSolve()
25 
26    Notes:
27    The available visualization contexts include
28 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
29 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
30          output where only the first processor opens
31          the file.  All other processors send their
32          data to the first processor to print.
33 
34    The user can open an alternative visualization context with
35    PetscViewerASCIIOpen() - output to a specified file.
36 
37    Level: beginner
38 
39 .keywords: SNES, view
40 
41 .seealso: PetscViewerASCIIOpen()
42 @*/
43 int SNESView(SNES snes,PetscViewer viewer)
44 {
45   SNES_KSP_EW_ConvCtx *kctx;
46   int                 ierr;
47   KSP                 ksp;
48   char                *type;
49   PetscTruth          isascii,isstring;
50 
51   PetscFunctionBegin;
52   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
53   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);
54   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
55   PetscCheckSameComm(snes,viewer);
56 
57   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
58   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
59   if (isascii) {
60     if (snes->prefix) {
61       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr);
62     } else {
63       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
64     }
65     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
66     if (type) {
67       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
68     } else {
69       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
70     }
71     if (snes->view) {
72       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
73       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
74       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
75     }
76     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
77     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",
78                  snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr);
79     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
80     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
81     if (snes->ksp_ewconv) {
82       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
83       if (kctx) {
84         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
85         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
86         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
87       }
88     }
89   } else if (isstring) {
90     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
91     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
92   }
93   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
94   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
95   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
96   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 /*
101   We retain a list of functions that also take SNES command
102   line options. These are called at the end SNESSetFromOptions()
103 */
104 #define MAXSETFROMOPTIONS 5
105 static int numberofsetfromoptions;
106 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "SNESAddOptionsChecker"
110 /*@C
111   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
112 
113   Not Collective
114 
115   Input Parameter:
116 . snescheck - function that checks for options
117 
118   Level: developer
119 
120 .seealso: SNESSetFromOptions()
121 @*/
122 int SNESAddOptionsChecker(int (*snescheck)(SNES))
123 {
124   PetscFunctionBegin;
125   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
126     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS);
127   }
128 
129   othersetfromoptions[numberofsetfromoptions++] = snescheck;
130   PetscFunctionReturn(0);
131 }
132 
133 #undef __FUNCT__
134 #define __FUNCT__ "SNESSetFromOptions"
135 /*@
136    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
137 
138    Collective on SNES
139 
140    Input Parameter:
141 .  snes - the SNES context
142 
143    Options Database Keys:
144 +  -snes_type <type> - ls, tr, umls, umtr, test
145 .  -snes_stol - convergence tolerance in terms of the norm
146                 of the change in the solution between steps
147 .  -snes_atol <atol> - absolute tolerance of residual norm
148 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
149 .  -snes_max_it <max_it> - maximum number of iterations
150 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
151 .  -snes_max_fail <max_fail> - maximum number of failures
152 .  -snes_trtol <trtol> - trust region tolerance
153 .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
154                                solver; hence iterations will continue until max_it
155                                or some other criterion is reached. Saves expense
156                                of convergence test
157 .  -snes_monitor - prints residual norm at each iteration
158 .  -snes_vecmonitor - plots solution at each iteration
159 .  -snes_vecmonitor_update - plots update to solution at each iteration
160 .  -snes_xmonitor - plots residual norm at each iteration
161 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
162 -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
163 
164     Options Database for Eisenstat-Walker method:
165 +  -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence
166 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
167 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
168 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
169 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
170 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
171 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
172 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
173 
174    Notes:
175    To see all options, run your program with the -help option or consult
176    the users manual.
177 
178    Level: beginner
179 
180 .keywords: SNES, nonlinear, set, options, database
181 
182 .seealso: SNESSetOptionsPrefix()
183 @*/
184 int SNESSetFromOptions(SNES snes)
185 {
186   KSP                 ksp;
187   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
188   PetscTruth          flg;
189   int                 ierr, i;
190   const char          *deft;
191   char                type[256];
192 
193   PetscFunctionBegin;
194   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
195 
196   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
197     if (snes->type_name) {
198       deft = snes->type_name;
199     } else {
200       deft = SNESLS;
201     }
202 
203     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
204     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
205     if (flg) {
206       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
207     } else if (!snes->type_name) {
208       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
209     }
210     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
211 
212     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
213     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
214 
215     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
216     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
217     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
218     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
219 
220     ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
221 
222     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
223     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
224     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
225     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
226     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
227     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
228     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
229 
230     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
231     if (flg) {snes->converged = 0;}
232     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
233     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
234     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
235     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
236     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
237     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
238     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
239     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
240     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
241     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
242     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
243     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
244     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
245     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
246     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
247     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
248 
249     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
250     if (flg) {
251       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
252       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
253     }
254 
255     for(i = 0; i < numberofsetfromoptions; i++) {
256       ierr = (*othersetfromoptions[i])(snes);                                                             CHKERRQ(ierr);
257     }
258 
259     if (snes->setfromoptions) {
260       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
261     }
262 
263   ierr = PetscOptionsEnd();CHKERRQ(ierr);
264 
265   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
266   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
267 
268   PetscFunctionReturn(0);
269 }
270 
271 
272 #undef __FUNCT__
273 #define __FUNCT__ "SNESSetApplicationContext"
274 /*@
275    SNESSetApplicationContext - Sets the optional user-defined context for
276    the nonlinear solvers.
277 
278    Collective on SNES
279 
280    Input Parameters:
281 +  snes - the SNES context
282 -  usrP - optional user context
283 
284    Level: intermediate
285 
286 .keywords: SNES, nonlinear, set, application, context
287 
288 .seealso: SNESGetApplicationContext()
289 @*/
290 int SNESSetApplicationContext(SNES snes,void *usrP)
291 {
292   PetscFunctionBegin;
293   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
294   snes->user		= usrP;
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "SNESGetApplicationContext"
300 /*@C
301    SNESGetApplicationContext - Gets the user-defined context for the
302    nonlinear solvers.
303 
304    Not Collective
305 
306    Input Parameter:
307 .  snes - SNES context
308 
309    Output Parameter:
310 .  usrP - user context
311 
312    Level: intermediate
313 
314 .keywords: SNES, nonlinear, get, application, context
315 
316 .seealso: SNESSetApplicationContext()
317 @*/
318 int SNESGetApplicationContext(SNES snes,void **usrP)
319 {
320   PetscFunctionBegin;
321   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
322   *usrP = snes->user;
323   PetscFunctionReturn(0);
324 }
325 
326 #undef __FUNCT__
327 #define __FUNCT__ "SNESGetIterationNumber"
328 /*@
329    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
330    at this time.
331 
332    Not Collective
333 
334    Input Parameter:
335 .  snes - SNES context
336 
337    Output Parameter:
338 .  iter - iteration number
339 
340    Notes:
341    For example, during the computation of iteration 2 this would return 1.
342 
343    This is useful for using lagged Jacobians (where one does not recompute the
344    Jacobian at each SNES iteration). For example, the code
345 .vb
346       ierr = SNESGetIterationNumber(snes,&it);
347       if (!(it % 2)) {
348         [compute Jacobian here]
349       }
350 .ve
351    can be used in your ComputeJacobian() function to cause the Jacobian to be
352    recomputed every second SNES iteration.
353 
354    Level: intermediate
355 
356 .keywords: SNES, nonlinear, get, iteration, number
357 @*/
358 int SNESGetIterationNumber(SNES snes,int* iter)
359 {
360   PetscFunctionBegin;
361   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
362   PetscValidIntPointer(iter,2);
363   *iter = snes->iter;
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "SNESGetFunctionNorm"
369 /*@
370    SNESGetFunctionNorm - Gets the norm of the current function that was set
371    with SNESSSetFunction().
372 
373    Collective on SNES
374 
375    Input Parameter:
376 .  snes - SNES context
377 
378    Output Parameter:
379 .  fnorm - 2-norm of function
380 
381    Level: intermediate
382 
383 .keywords: SNES, nonlinear, get, function, norm
384 
385 .seealso: SNESGetFunction()
386 @*/
387 int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
388 {
389   PetscFunctionBegin;
390   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
391   PetscValidScalarPointer(fnorm,2);
392   *fnorm = snes->norm;
393   PetscFunctionReturn(0);
394 }
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
398 /*@
399    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
400    attempted by the nonlinear solver.
401 
402    Not Collective
403 
404    Input Parameter:
405 .  snes - SNES context
406 
407    Output Parameter:
408 .  nfails - number of unsuccessful steps attempted
409 
410    Notes:
411    This counter is reset to zero for each successive call to SNESSolve().
412 
413    Level: intermediate
414 
415 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
416 @*/
417 int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
418 {
419   PetscFunctionBegin;
420   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
421   PetscValidIntPointer(nfails,2);
422   *nfails = snes->numFailures;
423   PetscFunctionReturn(0);
424 }
425 
426 #undef __FUNCT__
427 #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
428 /*@
429    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
430    attempted by the nonlinear solver before it gives up.
431 
432    Not Collective
433 
434    Input Parameters:
435 +  snes     - SNES context
436 -  maxFails - maximum of unsuccessful steps
437 
438    Level: intermediate
439 
440 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
441 @*/
442 int SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails)
443 {
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
446   snes->maxFailures = maxFails;
447   PetscFunctionReturn(0);
448 }
449 
450 #undef __FUNCT__
451 #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
452 /*@
453    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
454    attempted by the nonlinear solver before it gives up.
455 
456    Not Collective
457 
458    Input Parameter:
459 .  snes     - SNES context
460 
461    Output Parameter:
462 .  maxFails - maximum of unsuccessful steps
463 
464    Level: intermediate
465 
466 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
467 @*/
468 int SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails)
469 {
470   PetscFunctionBegin;
471   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
472   PetscValidIntPointer(maxFails,2);
473   *maxFails = snes->maxFailures;
474   PetscFunctionReturn(0);
475 }
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "SNESGetNumberLinearIterations"
479 /*@
480    SNESGetNumberLinearIterations - Gets the total number of linear iterations
481    used by the nonlinear solver.
482 
483    Not Collective
484 
485    Input Parameter:
486 .  snes - SNES context
487 
488    Output Parameter:
489 .  lits - number of linear iterations
490 
491    Notes:
492    This counter is reset to zero for each successive call to SNESSolve().
493 
494    Level: intermediate
495 
496 .keywords: SNES, nonlinear, get, number, linear, iterations
497 @*/
498 int SNESGetNumberLinearIterations(SNES snes,int* lits)
499 {
500   PetscFunctionBegin;
501   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
502   PetscValidIntPointer(lits,2);
503   *lits = snes->linear_its;
504   PetscFunctionReturn(0);
505 }
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "SNESGetKSP"
509 /*@C
510    SNESGetKSP - Returns the KSP context for a SNES solver.
511 
512    Not Collective, but if SNES object is parallel, then KSP object is parallel
513 
514    Input Parameter:
515 .  snes - the SNES context
516 
517    Output Parameter:
518 .  ksp - the KSP context
519 
520    Notes:
521    The user can then directly manipulate the KSP context to set various
522    options, etc.  Likewise, the user can then extract and manipulate the
523    KSP and PC contexts as well.
524 
525    Level: beginner
526 
527 .keywords: SNES, nonlinear, get, KSP, context
528 
529 .seealso: KSPGetPC()
530 @*/
531 int SNESGetKSP(SNES snes,KSP *ksp)
532 {
533   PetscFunctionBegin;
534   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
535   PetscValidPointer(ksp,2);
536   *ksp = snes->ksp;
537   PetscFunctionReturn(0);
538 }
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "SNESPublish_Petsc"
542 static int SNESPublish_Petsc(PetscObject obj)
543 {
544 #if defined(PETSC_HAVE_AMS)
545   SNES          v = (SNES) obj;
546   int          ierr;
547 #endif
548 
549   PetscFunctionBegin;
550 
551 #if defined(PETSC_HAVE_AMS)
552   /* if it is already published then return */
553   if (v->amem >=0) PetscFunctionReturn(0);
554 
555   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
556   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
557                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
558   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
559                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
560   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
561 #endif
562   PetscFunctionReturn(0);
563 }
564 
565 /* -----------------------------------------------------------*/
566 #undef __FUNCT__
567 #define __FUNCT__ "SNESCreate"
568 /*@C
569    SNESCreate - Creates a nonlinear solver context.
570 
571    Collective on MPI_Comm
572 
573    Input Parameters:
574 +  comm - MPI communicator
575 
576    Output Parameter:
577 .  outsnes - the new SNES context
578 
579    Options Database Keys:
580 +   -snes_mf - Activates default matrix-free Jacobian-vector products,
581                and no preconditioning matrix
582 .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
583                products, and a user-provided preconditioning matrix
584                as set by SNESSetJacobian()
585 -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
586 
587    Level: beginner
588 
589 .keywords: SNES, nonlinear, create, context
590 
591 .seealso: SNESSolve(), SNESDestroy(), SNES
592 @*/
593 int SNESCreate(MPI_Comm comm,SNES *outsnes)
594 {
595   int                 ierr;
596   SNES                snes;
597   SNES_KSP_EW_ConvCtx *kctx;
598 
599   PetscFunctionBegin;
600   PetscValidPointer(outsnes,1);
601   *outsnes = PETSC_NULL;
602 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
603   ierr = SNESInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
604 #endif
605 
606   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
607   PetscLogObjectCreate(snes);
608   snes->bops->publish     = SNESPublish_Petsc;
609   snes->max_its           = 50;
610   snes->max_funcs	  = 10000;
611   snes->norm		  = 0.0;
612   snes->rtol		  = 1.e-8;
613   snes->ttol              = 0.0;
614   snes->atol		  = 1.e-50;
615   snes->xtol		  = 1.e-8;
616   snes->deltatol	  = 1.e-12;
617   snes->nfuncs            = 0;
618   snes->numFailures       = 0;
619   snes->maxFailures       = 1;
620   snes->linear_its        = 0;
621   snes->numbermonitors    = 0;
622   snes->data              = 0;
623   snes->view              = 0;
624   snes->setupcalled       = 0;
625   snes->ksp_ewconv        = PETSC_FALSE;
626   snes->vwork             = 0;
627   snes->nwork             = 0;
628   snes->conv_hist_len     = 0;
629   snes->conv_hist_max     = 0;
630   snes->conv_hist         = PETSC_NULL;
631   snes->conv_hist_its     = PETSC_NULL;
632   snes->conv_hist_reset   = PETSC_TRUE;
633   snes->reason            = SNES_CONVERGED_ITERATING;
634 
635   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
636   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
637   PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
638   snes->kspconvctx  = (void*)kctx;
639   kctx->version     = 2;
640   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
641                              this was too large for some test cases */
642   kctx->rtol_last   = 0;
643   kctx->rtol_max    = .9;
644   kctx->gamma       = 1.0;
645   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
646   kctx->alpha       = kctx->alpha2;
647   kctx->threshold   = .1;
648   kctx->lresid_last = 0;
649   kctx->norm_last   = 0;
650 
651   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
652   PetscLogObjectParent(snes,snes->ksp)
653 
654   *outsnes = snes;
655   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
656   PetscFunctionReturn(0);
657 }
658 
659 /* --------------------------------------------------------------- */
660 #undef __FUNCT__
661 #define __FUNCT__ "SNESSetFunction"
662 /*@C
663    SNESSetFunction - Sets the function evaluation routine and function
664    vector for use by the SNES routines in solving systems of nonlinear
665    equations.
666 
667    Collective on SNES
668 
669    Input Parameters:
670 +  snes - the SNES context
671 .  func - function evaluation routine
672 .  r - vector to store function value
673 -  ctx - [optional] user-defined context for private data for the
674          function evaluation routine (may be PETSC_NULL)
675 
676    Calling sequence of func:
677 $    func (SNES snes,Vec x,Vec f,void *ctx);
678 
679 .  f - function vector
680 -  ctx - optional user-defined function context
681 
682    Notes:
683    The Newton-like methods typically solve linear systems of the form
684 $      f'(x) x = -f(x),
685    where f'(x) denotes the Jacobian matrix and f(x) is the function.
686 
687    Level: beginner
688 
689 .keywords: SNES, nonlinear, set, function
690 
691 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
692 @*/
693 int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
694 {
695   PetscFunctionBegin;
696   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
697   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
698   PetscCheckSameComm(snes,r);
699 
700   snes->computefunction     = func;
701   snes->vec_func            = snes->vec_func_always = r;
702   snes->funP                = ctx;
703   PetscFunctionReturn(0);
704 }
705 
706 #undef __FUNCT__
707 #define __FUNCT__ "SNESComputeFunction"
708 /*@
709    SNESComputeFunction - Calls the function that has been set with
710                          SNESSetFunction().
711 
712    Collective on SNES
713 
714    Input Parameters:
715 +  snes - the SNES context
716 -  x - input vector
717 
718    Output Parameter:
719 .  y - function vector, as set by SNESSetFunction()
720 
721    Notes:
722    SNESComputeFunction() is typically used within nonlinear solvers
723    implementations, so most users would not generally call this routine
724    themselves.
725 
726    Level: developer
727 
728 .keywords: SNES, nonlinear, compute, function
729 
730 .seealso: SNESSetFunction(), SNESGetFunction()
731 @*/
732 int SNESComputeFunction(SNES snes,Vec x,Vec y)
733 {
734   int    ierr;
735 
736   PetscFunctionBegin;
737   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
738   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
739   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
740   PetscCheckSameComm(snes,x);
741   PetscCheckSameComm(snes,y);
742 
743   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
744   PetscStackPush("SNES user function");
745   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
746   PetscStackPop;
747   snes->nfuncs++;
748   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "SNESComputeJacobian"
754 /*@
755    SNESComputeJacobian - Computes the Jacobian matrix that has been
756    set with SNESSetJacobian().
757 
758    Collective on SNES and Mat
759 
760    Input Parameters:
761 +  snes - the SNES context
762 -  x - input vector
763 
764    Output Parameters:
765 +  A - Jacobian matrix
766 .  B - optional preconditioning matrix
767 -  flag - flag indicating matrix structure
768 
769    Notes:
770    Most users should not need to explicitly call this routine, as it
771    is used internally within the nonlinear solvers.
772 
773    See KSPSetOperators() for important information about setting the
774    flag parameter.
775 
776    Level: developer
777 
778 .keywords: SNES, compute, Jacobian, matrix
779 
780 .seealso:  SNESSetJacobian(), KSPSetOperators()
781 @*/
782 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
783 {
784   int    ierr;
785 
786   PetscFunctionBegin;
787   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
788   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
789   PetscValidPointer(flg,5);
790   PetscCheckSameComm(snes,X);
791   if (!snes->computejacobian) PetscFunctionReturn(0);
792   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
793   *flg = DIFFERENT_NONZERO_PATTERN;
794   PetscStackPush("SNES user Jacobian function");
795   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
796   PetscStackPop;
797   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
798   /* make sure user returned a correct Jacobian and preconditioner */
799   PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
800   PetscValidHeaderSpecific(*B,MAT_COOKIE,4);
801   PetscFunctionReturn(0);
802 }
803 
804 #undef __FUNCT__
805 #define __FUNCT__ "SNESSetJacobian"
806 /*@C
807    SNESSetJacobian - Sets the function to compute Jacobian as well as the
808    location to store the matrix.
809 
810    Collective on SNES and Mat
811 
812    Input Parameters:
813 +  snes - the SNES context
814 .  A - Jacobian matrix
815 .  B - preconditioner matrix (usually same as the Jacobian)
816 .  func - Jacobian evaluation routine
817 -  ctx - [optional] user-defined context for private data for the
818          Jacobian evaluation routine (may be PETSC_NULL)
819 
820    Calling sequence of func:
821 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
822 
823 +  x - input vector
824 .  A - Jacobian matrix
825 .  B - preconditioner matrix, usually the same as A
826 .  flag - flag indicating information about the preconditioner matrix
827    structure (same as flag in KSPSetOperators())
828 -  ctx - [optional] user-defined Jacobian context
829 
830    Notes:
831    See KSPSetOperators() for important information about setting the flag
832    output parameter in the routine func().  Be sure to read this information!
833 
834    The routine func() takes Mat * as the matrix arguments rather than Mat.
835    This allows the Jacobian evaluation routine to replace A and/or B with a
836    completely new new matrix structure (not just different matrix elements)
837    when appropriate, for instance, if the nonzero structure is changing
838    throughout the global iterations.
839 
840    Level: beginner
841 
842 .keywords: SNES, nonlinear, set, Jacobian, matrix
843 
844 .seealso: KSPSetOperators(), SNESSetFunction()
845 @*/
846 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
847 {
848   int ierr;
849 
850   PetscFunctionBegin;
851   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
852   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
853   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
854   if (A) PetscCheckSameComm(snes,A);
855   if (B) PetscCheckSameComm(snes,B);
856   if (func) snes->computejacobian = func;
857   if (ctx)  snes->jacP            = ctx;
858   if (A) {
859     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
860     snes->jacobian = A;
861     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
862   }
863   if (B) {
864     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
865     snes->jacobian_pre = B;
866     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
867   }
868   PetscFunctionReturn(0);
869 }
870 
871 #undef __FUNCT__
872 #define __FUNCT__ "SNESGetJacobian"
873 /*@C
874    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
875    provided context for evaluating the Jacobian.
876 
877    Not Collective, but Mat object will be parallel if SNES object is
878 
879    Input Parameter:
880 .  snes - the nonlinear solver context
881 
882    Output Parameters:
883 +  A - location to stash Jacobian matrix (or PETSC_NULL)
884 .  B - location to stash preconditioner matrix (or PETSC_NULL)
885 .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
886 -  func - location to put Jacobian function (or PETSC_NULL)
887 
888    Level: advanced
889 
890 .seealso: SNESSetJacobian(), SNESComputeJacobian()
891 @*/
892 int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
893 {
894   PetscFunctionBegin;
895   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
896   if (A)    *A    = snes->jacobian;
897   if (B)    *B    = snes->jacobian_pre;
898   if (ctx)  *ctx  = snes->jacP;
899   if (func) *func = snes->computejacobian;
900   PetscFunctionReturn(0);
901 }
902 
903 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
904 extern int SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
905 
906 #undef __FUNCT__
907 #define __FUNCT__ "SNESSetUp"
908 /*@
909    SNESSetUp - Sets up the internal data structures for the later use
910    of a nonlinear solver.
911 
912    Collective on SNES
913 
914    Input Parameters:
915 +  snes - the SNES context
916 -  x - the solution vector
917 
918    Notes:
919    For basic use of the SNES solvers the user need not explicitly call
920    SNESSetUp(), since these actions will automatically occur during
921    the call to SNESSolve().  However, if one wishes to control this
922    phase separately, SNESSetUp() should be called after SNESCreate()
923    and optional routines of the form SNESSetXXX(), but before SNESSolve().
924 
925    Level: advanced
926 
927 .keywords: SNES, nonlinear, setup
928 
929 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
930 @*/
931 int SNESSetUp(SNES snes,Vec x)
932 {
933   int        ierr;
934   PetscTruth flg, iseqtr;
935 
936   PetscFunctionBegin;
937   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
938   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
939   PetscCheckSameComm(snes,x);
940   snes->vec_sol = snes->vec_sol_always = x;
941 
942   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
943   /*
944       This version replaces the user provided Jacobian matrix with a
945       matrix-free version but still employs the user-provided preconditioner matrix
946   */
947   if (flg) {
948     Mat J;
949     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
950     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
951     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n");
952     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
953     ierr = MatDestroy(J);CHKERRQ(ierr);
954   }
955 
956 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
957   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
958   if (flg) {
959     Mat J;
960     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
961     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
962     ierr = MatDestroy(J);CHKERRQ(ierr);
963   }
964 #endif
965 
966   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
967   /*
968       This version replaces both the user-provided Jacobian and the user-
969       provided preconditioner matrix with the default matrix free version.
970    */
971   if (flg) {
972     Mat  J;
973     KSP ksp;
974     PC   pc;
975 
976     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
977     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
978     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
979     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
980     ierr = MatDestroy(J);CHKERRQ(ierr);
981 
982     /* force no preconditioner */
983     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
984     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
985     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
986     if (!flg) {
987       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
988     }
989   }
990 
991   if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
992   if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
993   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
994   if (snes->vec_func == snes->vec_sol) {
995     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
996   }
997 
998   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
999   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
1000   if (snes->ksp_ewconv && !iseqtr) {
1001     KSP ksp;
1002     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1003     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1004   }
1005 
1006   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
1007   snes->setupcalled = 1;
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "SNESDestroy"
1013 /*@C
1014    SNESDestroy - Destroys the nonlinear solver context that was created
1015    with SNESCreate().
1016 
1017    Collective on SNES
1018 
1019    Input Parameter:
1020 .  snes - the SNES context
1021 
1022    Level: beginner
1023 
1024 .keywords: SNES, nonlinear, destroy
1025 
1026 .seealso: SNESCreate(), SNESSolve()
1027 @*/
1028 int SNESDestroy(SNES snes)
1029 {
1030   int i,ierr;
1031 
1032   PetscFunctionBegin;
1033   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1034   if (--snes->refct > 0) PetscFunctionReturn(0);
1035 
1036   /* if memory was published with AMS then destroy it */
1037   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1038 
1039   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1040   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
1041   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1042   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1043   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1044   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1045   for (i=0; i<snes->numbermonitors; i++) {
1046     if (snes->monitordestroy[i]) {
1047       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1048     }
1049   }
1050   PetscLogObjectDestroy((PetscObject)snes);
1051   PetscHeaderDestroy((PetscObject)snes);
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 /* ----------- Routines to set solver parameters ---------- */
1056 
1057 #undef __FUNCT__
1058 #define __FUNCT__ "SNESSetTolerances"
1059 /*@
1060    SNESSetTolerances - Sets various parameters used in convergence tests.
1061 
1062    Collective on SNES
1063 
1064    Input Parameters:
1065 +  snes - the SNES context
1066 .  atol - absolute convergence tolerance
1067 .  rtol - relative convergence tolerance
1068 .  stol -  convergence tolerance in terms of the norm
1069            of the change in the solution between steps
1070 .  maxit - maximum number of iterations
1071 -  maxf - maximum number of function evaluations
1072 
1073    Options Database Keys:
1074 +    -snes_atol <atol> - Sets atol
1075 .    -snes_rtol <rtol> - Sets rtol
1076 .    -snes_stol <stol> - Sets stol
1077 .    -snes_max_it <maxit> - Sets maxit
1078 -    -snes_max_funcs <maxf> - Sets maxf
1079 
1080    Notes:
1081    The default maximum number of iterations is 50.
1082    The default maximum number of function evaluations is 1000.
1083 
1084    Level: intermediate
1085 
1086 .keywords: SNES, nonlinear, set, convergence, tolerances
1087 
1088 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
1089 @*/
1090 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
1091 {
1092   PetscFunctionBegin;
1093   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1094   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1095   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1096   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1097   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1098   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1099   PetscFunctionReturn(0);
1100 }
1101 
1102 #undef __FUNCT__
1103 #define __FUNCT__ "SNESGetTolerances"
1104 /*@
1105    SNESGetTolerances - Gets various parameters used in convergence tests.
1106 
1107    Not Collective
1108 
1109    Input Parameters:
1110 +  snes - the SNES context
1111 .  atol - absolute convergence tolerance
1112 .  rtol - relative convergence tolerance
1113 .  stol -  convergence tolerance in terms of the norm
1114            of the change in the solution between steps
1115 .  maxit - maximum number of iterations
1116 -  maxf - maximum number of function evaluations
1117 
1118    Notes:
1119    The user can specify PETSC_NULL for any parameter that is not needed.
1120 
1121    Level: intermediate
1122 
1123 .keywords: SNES, nonlinear, get, convergence, tolerances
1124 
1125 .seealso: SNESSetTolerances()
1126 @*/
1127 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
1128 {
1129   PetscFunctionBegin;
1130   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1131   if (atol)  *atol  = snes->atol;
1132   if (rtol)  *rtol  = snes->rtol;
1133   if (stol)  *stol  = snes->xtol;
1134   if (maxit) *maxit = snes->max_its;
1135   if (maxf)  *maxf  = snes->max_funcs;
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 #undef __FUNCT__
1140 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1141 /*@
1142    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1143 
1144    Collective on SNES
1145 
1146    Input Parameters:
1147 +  snes - the SNES context
1148 -  tol - tolerance
1149 
1150    Options Database Key:
1151 .  -snes_trtol <tol> - Sets tol
1152 
1153    Level: intermediate
1154 
1155 .keywords: SNES, nonlinear, set, trust region, tolerance
1156 
1157 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
1158 @*/
1159 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1160 {
1161   PetscFunctionBegin;
1162   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1163   snes->deltatol = tol;
1164   PetscFunctionReturn(0);
1165 }
1166 
1167 /*
1168    Duplicate the lg monitors for SNES from KSP; for some reason with
1169    dynamic libraries things don't work under Sun4 if we just use
1170    macros instead of functions
1171 */
1172 #undef __FUNCT__
1173 #define __FUNCT__ "SNESLGMonitor"
1174 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1175 {
1176   int ierr;
1177 
1178   PetscFunctionBegin;
1179   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1180   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 #undef __FUNCT__
1185 #define __FUNCT__ "SNESLGMonitorCreate"
1186 int SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1187 {
1188   int ierr;
1189 
1190   PetscFunctionBegin;
1191   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1192   PetscFunctionReturn(0);
1193 }
1194 
1195 #undef __FUNCT__
1196 #define __FUNCT__ "SNESLGMonitorDestroy"
1197 int SNESLGMonitorDestroy(PetscDrawLG draw)
1198 {
1199   int ierr;
1200 
1201   PetscFunctionBegin;
1202   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 /* ------------ Routines to set performance monitoring options ----------- */
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "SNESSetMonitor"
1210 /*@C
1211    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
1212    iteration of the nonlinear solver to display the iteration's
1213    progress.
1214 
1215    Collective on SNES
1216 
1217    Input Parameters:
1218 +  snes - the SNES context
1219 .  func - monitoring routine
1220 .  mctx - [optional] user-defined context for private data for the
1221           monitor routine (use PETSC_NULL if no context is desitre)
1222 -  monitordestroy - [optional] routine that frees monitor context
1223           (may be PETSC_NULL)
1224 
1225    Calling sequence of func:
1226 $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1227 
1228 +    snes - the SNES context
1229 .    its - iteration number
1230 .    norm - 2-norm function value (may be estimated)
1231 -    mctx - [optional] monitoring context
1232 
1233    Options Database Keys:
1234 +    -snes_monitor        - sets SNESDefaultMonitor()
1235 .    -snes_xmonitor       - sets line graph monitor,
1236                             uses SNESLGMonitorCreate()
1237 _    -snes_cancelmonitors - cancels all monitors that have
1238                             been hardwired into a code by
1239                             calls to SNESSetMonitor(), but
1240                             does not cancel those set via
1241                             the options database.
1242 
1243    Notes:
1244    Several different monitoring routines may be set by calling
1245    SNESSetMonitor() multiple times; all will be called in the
1246    order in which they were set.
1247 
1248    Level: intermediate
1249 
1250 .keywords: SNES, nonlinear, set, monitor
1251 
1252 .seealso: SNESDefaultMonitor(), SNESClearMonitor()
1253 @*/
1254 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
1255 {
1256   PetscFunctionBegin;
1257   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1258   if (snes->numbermonitors >= MAXSNESMONITORS) {
1259     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1260   }
1261 
1262   snes->monitor[snes->numbermonitors]           = func;
1263   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1264   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1265   PetscFunctionReturn(0);
1266 }
1267 
1268 #undef __FUNCT__
1269 #define __FUNCT__ "SNESClearMonitor"
1270 /*@C
1271    SNESClearMonitor - Clears all the monitor functions for a SNES object.
1272 
1273    Collective on SNES
1274 
1275    Input Parameters:
1276 .  snes - the SNES context
1277 
1278    Options Database Key:
1279 .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1280     into a code by calls to SNESSetMonitor(), but does not cancel those
1281     set via the options database
1282 
1283    Notes:
1284    There is no way to clear one specific monitor from a SNES object.
1285 
1286    Level: intermediate
1287 
1288 .keywords: SNES, nonlinear, set, monitor
1289 
1290 .seealso: SNESDefaultMonitor(), SNESSetMonitor()
1291 @*/
1292 int SNESClearMonitor(SNES snes)
1293 {
1294   PetscFunctionBegin;
1295   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1296   snes->numbermonitors = 0;
1297   PetscFunctionReturn(0);
1298 }
1299 
1300 #undef __FUNCT__
1301 #define __FUNCT__ "SNESSetConvergenceTest"
1302 /*@C
1303    SNESSetConvergenceTest - Sets the function that is to be used
1304    to test for convergence of the nonlinear iterative solution.
1305 
1306    Collective on SNES
1307 
1308    Input Parameters:
1309 +  snes - the SNES context
1310 .  func - routine to test for convergence
1311 -  cctx - [optional] context for private data for the convergence routine
1312           (may be PETSC_NULL)
1313 
1314    Calling sequence of func:
1315 $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1316 
1317 +    snes - the SNES context
1318 .    cctx - [optional] convergence context
1319 .    reason - reason for convergence/divergence
1320 .    xnorm - 2-norm of current iterate
1321 .    gnorm - 2-norm of current step
1322 -    f - 2-norm of function
1323 
1324    Level: advanced
1325 
1326 .keywords: SNES, nonlinear, set, convergence, test
1327 
1328 .seealso: SNESConverged_LS(), SNESConverged_TR()
1329 @*/
1330 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
1331 {
1332   PetscFunctionBegin;
1333   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1334   (snes)->converged = func;
1335   (snes)->cnvP      = cctx;
1336   PetscFunctionReturn(0);
1337 }
1338 
1339 #undef __FUNCT__
1340 #define __FUNCT__ "SNESGetConvergedReason"
1341 /*@C
1342    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1343 
1344    Not Collective
1345 
1346    Input Parameter:
1347 .  snes - the SNES context
1348 
1349    Output Parameter:
1350 .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1351             manual pages for the individual convergence tests for complete lists
1352 
1353    Level: intermediate
1354 
1355    Notes: Can only be called after the call the SNESSolve() is complete.
1356 
1357 .keywords: SNES, nonlinear, set, convergence, test
1358 
1359 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1360 @*/
1361 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1362 {
1363   PetscFunctionBegin;
1364   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1365   PetscValidPointer(reason,2);
1366   *reason = snes->reason;
1367   PetscFunctionReturn(0);
1368 }
1369 
1370 #undef __FUNCT__
1371 #define __FUNCT__ "SNESSetConvergenceHistory"
1372 /*@
1373    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1374 
1375    Collective on SNES
1376 
1377    Input Parameters:
1378 +  snes - iterative context obtained from SNESCreate()
1379 .  a   - array to hold history
1380 .  its - integer array holds the number of linear iterations for each solve.
1381 .  na  - size of a and its
1382 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1383            else it continues storing new values for new nonlinear solves after the old ones
1384 
1385    Notes:
1386    If set, this array will contain the function norms computed
1387    at each step.
1388 
1389    This routine is useful, e.g., when running a code for purposes
1390    of accurate performance monitoring, when no I/O should be done
1391    during the section of code that is being timed.
1392 
1393    Level: intermediate
1394 
1395 .keywords: SNES, set, convergence, history
1396 
1397 .seealso: SNESGetConvergenceHistory()
1398 
1399 @*/
1400 int SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset)
1401 {
1402   PetscFunctionBegin;
1403   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1404   if (na) PetscValidScalarPointer(a,2);
1405   snes->conv_hist       = a;
1406   snes->conv_hist_its   = its;
1407   snes->conv_hist_max   = na;
1408   snes->conv_hist_reset = reset;
1409   PetscFunctionReturn(0);
1410 }
1411 
1412 #undef __FUNCT__
1413 #define __FUNCT__ "SNESGetConvergenceHistory"
1414 /*@C
1415    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1416 
1417    Collective on SNES
1418 
1419    Input Parameter:
1420 .  snes - iterative context obtained from SNESCreate()
1421 
1422    Output Parameters:
1423 .  a   - array to hold history
1424 .  its - integer array holds the number of linear iterations (or
1425          negative if not converged) for each solve.
1426 -  na  - size of a and its
1427 
1428    Notes:
1429     The calling sequence for this routine in Fortran is
1430 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1431 
1432    This routine is useful, e.g., when running a code for purposes
1433    of accurate performance monitoring, when no I/O should be done
1434    during the section of code that is being timed.
1435 
1436    Level: intermediate
1437 
1438 .keywords: SNES, get, convergence, history
1439 
1440 .seealso: SNESSetConvergencHistory()
1441 
1442 @*/
1443 int SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na)
1444 {
1445   PetscFunctionBegin;
1446   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1447   if (a)   *a   = snes->conv_hist;
1448   if (its) *its = snes->conv_hist_its;
1449   if (na) *na   = snes->conv_hist_len;
1450   PetscFunctionReturn(0);
1451 }
1452 
1453 #undef __FUNCT__
1454 #define __FUNCT__ "SNESSetRhsBC"
1455 /*@C
1456   SNESSetRhsBC - Sets the function which applies boundary conditions
1457   to the Rhs of each system.
1458 
1459   Collective on SNES
1460 
1461   Input Parameters:
1462 . snes - The nonlinear solver context
1463 . func - The function
1464 
1465   Calling sequence of func:
1466 . func (SNES snes, Vec rhs, void *ctx);
1467 
1468 . rhs - The current rhs vector
1469 . ctx - The user-context
1470 
1471   Level: intermediate
1472 
1473 .keywords: SNES, Rhs, boundary conditions
1474 .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
1475 @*/
1476 int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *))
1477 {
1478   PetscFunctionBegin;
1479   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
1480   snes->applyrhsbc = func;
1481   PetscFunctionReturn(0);
1482 }
1483 
1484 #undef __FUNCT__
1485 #define __FUNCT__ "SNESDefaultRhsBC"
1486 /*@
1487   SNESDefaultRhsBC - The default boundary condition function which does nothing.
1488 
1489   Not collective
1490 
1491   Input Parameters:
1492 . snes - The nonlinear solver context
1493 . rhs  - The Rhs
1494 . ctx  - The user-context
1495 
1496   Level: beginner
1497 
1498 .keywords: SNES, Rhs, boundary conditions
1499 .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
1500 @*/
1501 int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
1502 {
1503   PetscFunctionBegin;
1504   PetscFunctionReturn(0);
1505 }
1506 
1507 #undef __FUNCT__
1508 #define __FUNCT__ "SNESSetSolutionBC"
1509 /*@C
1510   SNESSetSolutionBC - Sets the function which applies boundary conditions
1511   to the solution of each system.
1512 
1513   Collective on SNES
1514 
1515   Input Parameters:
1516 . snes - The nonlinear solver context
1517 . func - The function
1518 
1519   Calling sequence of func:
1520 . func (SNES snes, Vec rsol, void *ctx);
1521 
1522 . sol - The current solution vector
1523 . ctx - The user-context
1524 
1525   Level: intermediate
1526 
1527 .keywords: SNES, solution, boundary conditions
1528 .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
1529 @*/
1530 int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *))
1531 {
1532   PetscFunctionBegin;
1533   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
1534   snes->applysolbc = func;
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 #undef __FUNCT__
1539 #define __FUNCT__ "SNESDefaultSolutionBC"
1540 /*@
1541   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
1542 
1543   Not collective
1544 
1545   Input Parameters:
1546 . snes - The nonlinear solver context
1547 . sol  - The solution
1548 . ctx  - The user-context
1549 
1550   Level: beginner
1551 
1552 .keywords: SNES, solution, boundary conditions
1553 .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
1554 @*/
1555 int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
1556 {
1557   PetscFunctionBegin;
1558   PetscFunctionReturn(0);
1559 }
1560 
1561 #undef __FUNCT__
1562 #define __FUNCT__ "SNESSetUpdate"
1563 /*@C
1564   SNESSetUpdate - Sets the general-purpose update function called
1565   at the beginning of every step of the iteration.
1566 
1567   Collective on SNES
1568 
1569   Input Parameters:
1570 . snes - The nonlinear solver context
1571 . func - The function
1572 
1573   Calling sequence of func:
1574 . func (TS ts, int step);
1575 
1576 . step - The current step of the iteration
1577 
1578   Level: intermediate
1579 
1580 .keywords: SNES, update
1581 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
1582 @*/
1583 int SNESSetUpdate(SNES snes, int (*func)(SNES, int))
1584 {
1585   PetscFunctionBegin;
1586   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
1587   snes->update = func;
1588   PetscFunctionReturn(0);
1589 }
1590 
1591 #undef __FUNCT__
1592 #define __FUNCT__ "SNESDefaultUpdate"
1593 /*@
1594   SNESDefaultUpdate - The default update function which does nothing.
1595 
1596   Not collective
1597 
1598   Input Parameters:
1599 . snes - The nonlinear solver context
1600 . step - The current step of the iteration
1601 
1602   Level: intermediate
1603 
1604 .keywords: SNES, update
1605 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
1606 @*/
1607 int SNESDefaultUpdate(SNES snes, int step)
1608 {
1609   PetscFunctionBegin;
1610   PetscFunctionReturn(0);
1611 }
1612 
1613 #undef __FUNCT__
1614 #define __FUNCT__ "SNESScaleStep_Private"
1615 /*
1616    SNESScaleStep_Private - Scales a step so that its length is less than the
1617    positive parameter delta.
1618 
1619     Input Parameters:
1620 +   snes - the SNES context
1621 .   y - approximate solution of linear system
1622 .   fnorm - 2-norm of current function
1623 -   delta - trust region size
1624 
1625     Output Parameters:
1626 +   gpnorm - predicted function norm at the new point, assuming local
1627     linearization.  The value is zero if the step lies within the trust
1628     region, and exceeds zero otherwise.
1629 -   ynorm - 2-norm of the step
1630 
1631     Note:
1632     For non-trust region methods such as SNESLS, the parameter delta
1633     is set to be the maximum allowable step size.
1634 
1635 .keywords: SNES, nonlinear, scale, step
1636 */
1637 int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
1638 {
1639   PetscReal   nrm;
1640   PetscScalar cnorm;
1641   int         ierr;
1642 
1643   PetscFunctionBegin;
1644   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1645   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1646   PetscCheckSameComm(snes,y);
1647 
1648   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1649   if (nrm > *delta) {
1650      nrm = *delta/nrm;
1651      *gpnorm = (1.0 - nrm)*(*fnorm);
1652      cnorm = nrm;
1653      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
1654      *ynorm = *delta;
1655   } else {
1656      *gpnorm = 0.0;
1657      *ynorm = nrm;
1658   }
1659   PetscFunctionReturn(0);
1660 }
1661 
1662 #undef __FUNCT__
1663 #define __FUNCT__ "SNESSolve"
1664 /*@
1665    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
1666    SNESCreate() and optional routines of the form SNESSetXXX().
1667 
1668    Collective on SNES
1669 
1670    Input Parameters:
1671 +  snes - the SNES context
1672 -  x - the solution vector
1673 
1674    Output Parameter:
1675 .  its - number of iterations until termination
1676 
1677    Notes:
1678    The user should initialize the vector,x, with the initial guess
1679    for the nonlinear solve prior to calling SNESSolve.  In particular,
1680    to employ an initial guess of zero, the user should explicitly set
1681    this vector to zero by calling VecSet().
1682 
1683    Level: beginner
1684 
1685 .keywords: SNES, nonlinear, solve
1686 
1687 .seealso: SNESCreate(), SNESDestroy()
1688 @*/
1689 int SNESSolve(SNES snes,Vec x,int *its)
1690 {
1691   int        ierr;
1692   PetscTruth flg;
1693 
1694   PetscFunctionBegin;
1695   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1696   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1697   PetscCheckSameComm(snes,x);
1698   PetscValidIntPointer(its,3);
1699   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
1700 
1701   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
1702   else {snes->vec_sol = snes->vec_sol_always = x;}
1703   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
1704   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1705   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1706   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
1707   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1708   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
1709   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); }
1710   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1711   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
1712   PetscFunctionReturn(0);
1713 }
1714 
1715 /* --------- Internal routines for SNES Package --------- */
1716 
1717 #undef __FUNCT__
1718 #define __FUNCT__ "SNESSetType"
1719 /*@C
1720    SNESSetType - Sets the method for the nonlinear solver.
1721 
1722    Collective on SNES
1723 
1724    Input Parameters:
1725 +  snes - the SNES context
1726 -  type - a known method
1727 
1728    Options Database Key:
1729 .  -snes_type <type> - Sets the method; use -help for a list
1730    of available methods (for instance, ls or tr)
1731 
1732    Notes:
1733    See "petsc/include/petscsnes.h" for available methods (for instance)
1734 +    SNESLS - Newton's method with line search
1735      (systems of nonlinear equations)
1736 .    SNESTR - Newton's method with trust region
1737      (systems of nonlinear equations)
1738 
1739   Normally, it is best to use the SNESSetFromOptions() command and then
1740   set the SNES solver type from the options database rather than by using
1741   this routine.  Using the options database provides the user with
1742   maximum flexibility in evaluating the many nonlinear solvers.
1743   The SNESSetType() routine is provided for those situations where it
1744   is necessary to set the nonlinear solver independently of the command
1745   line or options database.  This might be the case, for example, when
1746   the choice of solver changes during the execution of the program,
1747   and the user's application is taking responsibility for choosing the
1748   appropriate method.
1749 
1750   Level: intermediate
1751 
1752 .keywords: SNES, set, type
1753 
1754 .seealso: SNESType, SNESCreate()
1755 
1756 @*/
1757 int SNESSetType(SNES snes,const SNESType type)
1758 {
1759   int        ierr,(*r)(SNES);
1760   PetscTruth match;
1761 
1762   PetscFunctionBegin;
1763   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1764   PetscValidCharPointer(type,2);
1765 
1766   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
1767   if (match) PetscFunctionReturn(0);
1768 
1769   if (snes->setupcalled) {
1770     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
1771     snes->data = 0;
1772   }
1773 
1774   /* Get the function pointers for the iterative method requested */
1775   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1776 
1777   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1778 
1779   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
1780 
1781   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
1782   snes->data = 0;
1783   ierr = (*r)(snes);CHKERRQ(ierr);
1784   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
1785 
1786   PetscFunctionReturn(0);
1787 }
1788 
1789 
1790 /* --------------------------------------------------------------------- */
1791 #undef __FUNCT__
1792 #define __FUNCT__ "SNESRegisterDestroy"
1793 /*@C
1794    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1795    registered by SNESRegisterDynamic().
1796 
1797    Not Collective
1798 
1799    Level: advanced
1800 
1801 .keywords: SNES, nonlinear, register, destroy
1802 
1803 .seealso: SNESRegisterAll(), SNESRegisterAll()
1804 @*/
1805 int SNESRegisterDestroy(void)
1806 {
1807   int ierr;
1808 
1809   PetscFunctionBegin;
1810   if (SNESList) {
1811     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
1812     SNESList = 0;
1813   }
1814   SNESRegisterAllCalled = PETSC_FALSE;
1815   PetscFunctionReturn(0);
1816 }
1817 
1818 #undef __FUNCT__
1819 #define __FUNCT__ "SNESGetType"
1820 /*@C
1821    SNESGetType - Gets the SNES method type and name (as a string).
1822 
1823    Not Collective
1824 
1825    Input Parameter:
1826 .  snes - nonlinear solver context
1827 
1828    Output Parameter:
1829 .  type - SNES method (a character string)
1830 
1831    Level: intermediate
1832 
1833 .keywords: SNES, nonlinear, get, type, name
1834 @*/
1835 int SNESGetType(SNES snes,SNESType *type)
1836 {
1837   PetscFunctionBegin;
1838   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1839   PetscValidPointer(type,2);
1840   *type = snes->type_name;
1841   PetscFunctionReturn(0);
1842 }
1843 
1844 #undef __FUNCT__
1845 #define __FUNCT__ "SNESGetSolution"
1846 /*@C
1847    SNESGetSolution - Returns the vector where the approximate solution is
1848    stored.
1849 
1850    Not Collective, but Vec is parallel if SNES is parallel
1851 
1852    Input Parameter:
1853 .  snes - the SNES context
1854 
1855    Output Parameter:
1856 .  x - the solution
1857 
1858    Level: advanced
1859 
1860 .keywords: SNES, nonlinear, get, solution
1861 
1862 .seealso: SNESGetFunction(), SNESGetSolutionUpdate()
1863 @*/
1864 int SNESGetSolution(SNES snes,Vec *x)
1865 {
1866   PetscFunctionBegin;
1867   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1868   PetscValidPointer(x,2);
1869   *x = snes->vec_sol_always;
1870   PetscFunctionReturn(0);
1871 }
1872 
1873 #undef __FUNCT__
1874 #define __FUNCT__ "SNESGetSolutionUpdate"
1875 /*@C
1876    SNESGetSolutionUpdate - Returns the vector where the solution update is
1877    stored.
1878 
1879    Not Collective, but Vec is parallel if SNES is parallel
1880 
1881    Input Parameter:
1882 .  snes - the SNES context
1883 
1884    Output Parameter:
1885 .  x - the solution update
1886 
1887    Level: advanced
1888 
1889 .keywords: SNES, nonlinear, get, solution, update
1890 
1891 .seealso: SNESGetSolution(), SNESGetFunction
1892 @*/
1893 int SNESGetSolutionUpdate(SNES snes,Vec *x)
1894 {
1895   PetscFunctionBegin;
1896   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1897   PetscValidPointer(x,2);
1898   *x = snes->vec_sol_update_always;
1899   PetscFunctionReturn(0);
1900 }
1901 
1902 #undef __FUNCT__
1903 #define __FUNCT__ "SNESGetFunction"
1904 /*@C
1905    SNESGetFunction - Returns the vector where the function is stored.
1906 
1907    Not Collective, but Vec is parallel if SNES is parallel
1908 
1909    Input Parameter:
1910 .  snes - the SNES context
1911 
1912    Output Parameter:
1913 +  r - the function (or PETSC_NULL)
1914 .  ctx - the function context (or PETSC_NULL)
1915 -  func - the function (or PETSC_NULL)
1916 
1917    Level: advanced
1918 
1919 .keywords: SNES, nonlinear, get, function
1920 
1921 .seealso: SNESSetFunction(), SNESGetSolution()
1922 @*/
1923 int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
1924 {
1925   PetscFunctionBegin;
1926   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1927   if (r)    *r    = snes->vec_func_always;
1928   if (ctx)  *ctx  = snes->funP;
1929   if (func) *func = snes->computefunction;
1930   PetscFunctionReturn(0);
1931 }
1932 
1933 #undef __FUNCT__
1934 #define __FUNCT__ "SNESSetOptionsPrefix"
1935 /*@C
1936    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1937    SNES options in the database.
1938 
1939    Collective on SNES
1940 
1941    Input Parameter:
1942 +  snes - the SNES context
1943 -  prefix - the prefix to prepend to all option names
1944 
1945    Notes:
1946    A hyphen (-) must NOT be given at the beginning of the prefix name.
1947    The first character of all runtime options is AUTOMATICALLY the hyphen.
1948 
1949    Level: advanced
1950 
1951 .keywords: SNES, set, options, prefix, database
1952 
1953 .seealso: SNESSetFromOptions()
1954 @*/
1955 int SNESSetOptionsPrefix(SNES snes,const char prefix[])
1956 {
1957   int ierr;
1958 
1959   PetscFunctionBegin;
1960   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1961   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
1962   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
1963   PetscFunctionReturn(0);
1964 }
1965 
1966 #undef __FUNCT__
1967 #define __FUNCT__ "SNESAppendOptionsPrefix"
1968 /*@C
1969    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1970    SNES options in the database.
1971 
1972    Collective on SNES
1973 
1974    Input Parameters:
1975 +  snes - the SNES context
1976 -  prefix - the prefix to prepend to all option names
1977 
1978    Notes:
1979    A hyphen (-) must NOT be given at the beginning of the prefix name.
1980    The first character of all runtime options is AUTOMATICALLY the hyphen.
1981 
1982    Level: advanced
1983 
1984 .keywords: SNES, append, options, prefix, database
1985 
1986 .seealso: SNESGetOptionsPrefix()
1987 @*/
1988 int SNESAppendOptionsPrefix(SNES snes,const char prefix[])
1989 {
1990   int ierr;
1991 
1992   PetscFunctionBegin;
1993   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1994   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
1995   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
1996   PetscFunctionReturn(0);
1997 }
1998 
1999 #undef __FUNCT__
2000 #define __FUNCT__ "SNESGetOptionsPrefix"
2001 /*@C
2002    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2003    SNES options in the database.
2004 
2005    Not Collective
2006 
2007    Input Parameter:
2008 .  snes - the SNES context
2009 
2010    Output Parameter:
2011 .  prefix - pointer to the prefix string used
2012 
2013    Notes: On the fortran side, the user should pass in a string 'prifix' of
2014    sufficient length to hold the prefix.
2015 
2016    Level: advanced
2017 
2018 .keywords: SNES, get, options, prefix, database
2019 
2020 .seealso: SNESAppendOptionsPrefix()
2021 @*/
2022 int SNESGetOptionsPrefix(SNES snes,char *prefix[])
2023 {
2024   int ierr;
2025 
2026   PetscFunctionBegin;
2027   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2028   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2029   PetscFunctionReturn(0);
2030 }
2031 
2032 
2033 #undef __FUNCT__
2034 #define __FUNCT__ "SNESRegister"
2035 /*@C
2036   SNESRegister - See SNESRegisterDynamic()
2037 
2038   Level: advanced
2039 @*/
2040 int SNESRegister(const char sname[],const char path[],const char name[],int (*function)(SNES))
2041 {
2042   char fullname[256];
2043   int  ierr;
2044 
2045   PetscFunctionBegin;
2046   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2047   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2048   PetscFunctionReturn(0);
2049 }
2050 
2051 #undef __FUNCT__
2052 #define __FUNCT__ "SNESTestLocalMin"
2053 int SNESTestLocalMin(SNES snes)
2054 {
2055   int         ierr,N,i,j;
2056   Vec         u,uh,fh;
2057   PetscScalar value;
2058   PetscReal   norm;
2059 
2060   PetscFunctionBegin;
2061   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2062   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2063   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2064 
2065   /* currently only works for sequential */
2066   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2067   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2068   for (i=0; i<N; i++) {
2069     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2070     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr);
2071     for (j=-10; j<11; j++) {
2072       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2073       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2074       ierr  = (*snes->computefunction)(snes,uh,fh,snes->funP);CHKERRQ(ierr);
2075       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2076       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %d %18.16e\n",j,norm);CHKERRQ(ierr);
2077       value = -value;
2078       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2079     }
2080   }
2081   ierr = VecDestroy(uh);CHKERRQ(ierr);
2082   ierr = VecDestroy(fh);CHKERRQ(ierr);
2083   PetscFunctionReturn(0);
2084 }
2085