xref: /petsc/src/snes/interface/snes.c (revision d32f9abdbc052d6e1fd06679b17a55415c3aae30)
1 #define PETSCSNES_DLL
2 
3 #include "include/private/snesimpl.h"      /*I "petscsnes.h"  I*/
4 
5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6 PetscFList SNESList              = PETSC_NULL;
7 
8 /* Logging support */
9 PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE;
10 PetscLogEvent  SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval;
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "SNESSetFunctionDomainError"
14 /*@
15    SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not
16      in the functions domain. For example, negative pressure.
17 
18    Collective on SNES
19 
20    Input Parameters:
21 .  SNES - the SNES context
22 
23    Level: advanced
24 
25 .keywords: SNES, view
26 
27 .seealso: SNESCreate(), SNESSetFunction()
28 @*/
29 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionDomainError(SNES snes)
30 {
31   PetscFunctionBegin;
32   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
33   snes->domainerror = PETSC_TRUE;
34   PetscFunctionReturn(0);
35 }
36 
37 #undef __FUNCT__
38 #define __FUNCT__ "SNESView"
39 /*@C
40    SNESView - Prints the SNES data structure.
41 
42    Collective on SNES
43 
44    Input Parameters:
45 +  SNES - the SNES context
46 -  viewer - visualization context
47 
48    Options Database Key:
49 .  -snes_view - Calls SNESView() at end of SNESSolve()
50 
51    Notes:
52    The available visualization contexts include
53 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
54 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
55          output where only the first processor opens
56          the file.  All other processors send their
57          data to the first processor to print.
58 
59    The user can open an alternative visualization context with
60    PetscViewerASCIIOpen() - output to a specified file.
61 
62    Level: beginner
63 
64 .keywords: SNES, view
65 
66 .seealso: PetscViewerASCIIOpen()
67 @*/
68 PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer)
69 {
70   SNESKSPEW           *kctx;
71   PetscErrorCode      ierr;
72   KSP                 ksp;
73   const SNESType      type;
74   PetscTruth          iascii,isstring;
75 
76   PetscFunctionBegin;
77   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
78   if (!viewer) {
79     ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
80   }
81   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
82   PetscCheckSameComm(snes,1,viewer,2);
83 
84   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
85   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
86   if (iascii) {
87     if (((PetscObject)snes)->prefix) {
88       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr);
89     } else {
90       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
91     }
92     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
93     if (type) {
94       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
95     } else {
96       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
97     }
98     if (snes->ops->view) {
99       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
100       ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);
101       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
102     }
103     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
104     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
105                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
106     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
107     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
108     if (snes->ksp_ewconv) {
109       kctx = (SNESKSPEW *)snes->kspconvctx;
110       if (kctx) {
111         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
112         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
113         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
114       }
115     }
116   } else if (isstring) {
117     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
118     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
119   }
120   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
121   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
122   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
123   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
124   PetscFunctionReturn(0);
125 }
126 
127 /*
128   We retain a list of functions that also take SNES command
129   line options. These are called at the end SNESSetFromOptions()
130 */
131 #define MAXSETFROMOPTIONS 5
132 static PetscInt numberofsetfromoptions;
133 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "SNESAddOptionsChecker"
137 /*@C
138   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
139 
140   Not Collective
141 
142   Input Parameter:
143 . snescheck - function that checks for options
144 
145   Level: developer
146 
147 .seealso: SNESSetFromOptions()
148 @*/
149 PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
150 {
151   PetscFunctionBegin;
152   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
153     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
154   }
155   othersetfromoptions[numberofsetfromoptions++] = snescheck;
156   PetscFunctionReturn(0);
157 }
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "SNESSetFromOptions"
161 /*@
162    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
163 
164    Collective on SNES
165 
166    Input Parameter:
167 .  snes - the SNES context
168 
169    Options Database Keys:
170 +  -snes_type <type> - ls, tr, umls, umtr, test
171 .  -snes_stol - convergence tolerance in terms of the norm
172                 of the change in the solution between steps
173 .  -snes_atol <abstol> - absolute tolerance of residual norm
174 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
175 .  -snes_max_it <max_it> - maximum number of iterations
176 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
177 .  -snes_max_fail <max_fail> - maximum number of failures
178 .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
179 .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
180 .  -snes_trtol <trtol> - trust region tolerance
181 .  -snes_no_convergence_test - skip convergence test in nonlinear
182                                solver; hence iterations will continue until max_it
183                                or some other criterion is reached. Saves expense
184                                of convergence test
185 .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
186                                        filename given prints to stdout
187 .  -snes_monitor_solution - plots solution at each iteration
188 .  -snes_monitor_residual - plots residual (not its norm) at each iteration
189 .  -snes_monitor_solution_update - plots update to solution at each iteration
190 .  -snes_monitor_draw - plots residual norm at each iteration
191 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
192 .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
193 -  -snes_converged_reason - print the reason for convergence/divergence after each solve
194 
195     Options Database for Eisenstat-Walker method:
196 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
197 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
198 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
199 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
200 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
201 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
202 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
203 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
204 
205    Notes:
206    To see all options, run your program with the -help option or consult
207    the users manual.
208 
209    Level: beginner
210 
211 .keywords: SNES, nonlinear, set, options, database
212 
213 .seealso: SNESSetOptionsPrefix()
214 @*/
215 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes)
216 {
217   PetscTruth              flg;
218   PetscInt                i,indx,lag;
219   const char              *deft = SNESLS;
220   const char              *convtests[] = {"default","skip"};
221   SNESKSPEW               *kctx = NULL;
222   char                    type[256], monfilename[PETSC_MAX_PATH_LEN];
223   PetscViewerASCIIMonitor monviewer;
224   PetscErrorCode          ierr;
225 
226   PetscFunctionBegin;
227   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
228 
229   ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
230     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
231     if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; }
232     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
233     if (flg) {
234       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
235     } else if (!((PetscObject)snes)->type_name) {
236       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
237     }
238     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
239 
240     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
241     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
242 
243     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
244     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
245     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
246     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
247     ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr);
248 
249     ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr);
250     if (flg) {
251       ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr);
252     }
253 
254     ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr);
255     if (flg) {
256       switch (indx) {
257       case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break;
258       case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);    break;
259       }
260     }
261 
262     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
263     if (flg) {
264       snes->printreason = PETSC_TRUE;
265     }
266 
267     kctx = (SNESKSPEW *)snes->kspconvctx;
268 
269     ierr = PetscOptionsTruth("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr);
270 
271     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
272     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
273     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
274     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
275     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
276     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
277     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
278 
279     ierr = PetscOptionsName("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",&flg);CHKERRQ(ierr);
280     if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);}
281 
282     ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
283     if (flg) {
284       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
285       ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
286     }
287 
288     ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
289     if (flg) {
290       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
291       ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr);
292     }
293 
294     ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
295     if (flg) {
296       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
297       ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
298     }
299 
300     ierr = PetscOptionsName("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",&flg);CHKERRQ(ierr);
301     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);}
302     ierr = PetscOptionsName("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",&flg);CHKERRQ(ierr);
303     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);}
304     ierr = PetscOptionsName("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",&flg);CHKERRQ(ierr);
305     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);}
306     ierr = PetscOptionsName("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",&flg);CHKERRQ(ierr);
307     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
308 
309     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
310     if (flg) {
311       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
312       ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr);
313     }
314 
315     for(i = 0; i < numberofsetfromoptions; i++) {
316       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
317     }
318 
319     if (snes->ops->setfromoptions) {
320       ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr);
321     }
322   ierr = PetscOptionsEnd();CHKERRQ(ierr);
323 
324   ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr);
325 
326   PetscFunctionReturn(0);
327 }
328 
329 
330 #undef __FUNCT__
331 #define __FUNCT__ "SNESSetApplicationContext"
332 /*@
333    SNESSetApplicationContext - Sets the optional user-defined context for
334    the nonlinear solvers.
335 
336    Collective on SNES
337 
338    Input Parameters:
339 +  snes - the SNES context
340 -  usrP - optional user context
341 
342    Level: intermediate
343 
344 .keywords: SNES, nonlinear, set, application, context
345 
346 .seealso: SNESGetApplicationContext()
347 @*/
348 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP)
349 {
350   PetscFunctionBegin;
351   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
352   snes->user		= usrP;
353   PetscFunctionReturn(0);
354 }
355 
356 #undef __FUNCT__
357 #define __FUNCT__ "SNESGetApplicationContext"
358 /*@C
359    SNESGetApplicationContext - Gets the user-defined context for the
360    nonlinear solvers.
361 
362    Not Collective
363 
364    Input Parameter:
365 .  snes - SNES context
366 
367    Output Parameter:
368 .  usrP - user context
369 
370    Level: intermediate
371 
372 .keywords: SNES, nonlinear, get, application, context
373 
374 .seealso: SNESSetApplicationContext()
375 @*/
376 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP)
377 {
378   PetscFunctionBegin;
379   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
380   *usrP = snes->user;
381   PetscFunctionReturn(0);
382 }
383 
384 #undef __FUNCT__
385 #define __FUNCT__ "SNESGetIterationNumber"
386 /*@
387    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
388    at this time.
389 
390    Not Collective
391 
392    Input Parameter:
393 .  snes - SNES context
394 
395    Output Parameter:
396 .  iter - iteration number
397 
398    Notes:
399    For example, during the computation of iteration 2 this would return 1.
400 
401    This is useful for using lagged Jacobians (where one does not recompute the
402    Jacobian at each SNES iteration). For example, the code
403 .vb
404       ierr = SNESGetIterationNumber(snes,&it);
405       if (!(it % 2)) {
406         [compute Jacobian here]
407       }
408 .ve
409    can be used in your ComputeJacobian() function to cause the Jacobian to be
410    recomputed every second SNES iteration.
411 
412    Level: intermediate
413 
414 .keywords: SNES, nonlinear, get, iteration, number,
415 
416 .seealso:   SNESGetFunctionNorm(), SNESGetLinearSolveIterations()
417 @*/
418 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter)
419 {
420   PetscFunctionBegin;
421   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
422   PetscValidIntPointer(iter,2);
423   *iter = snes->iter;
424   PetscFunctionReturn(0);
425 }
426 
427 #undef __FUNCT__
428 #define __FUNCT__ "SNESGetFunctionNorm"
429 /*@
430    SNESGetFunctionNorm - Gets the norm of the current function that was set
431    with SNESSSetFunction().
432 
433    Collective on SNES
434 
435    Input Parameter:
436 .  snes - SNES context
437 
438    Output Parameter:
439 .  fnorm - 2-norm of function
440 
441    Level: intermediate
442 
443 .keywords: SNES, nonlinear, get, function, norm
444 
445 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations()
446 @*/
447 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscReal *fnorm)
448 {
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
451   PetscValidScalarPointer(fnorm,2);
452   *fnorm = snes->norm;
453   PetscFunctionReturn(0);
454 }
455 
456 #undef __FUNCT__
457 #define __FUNCT__ "SNESGetNonlinearStepFailures"
458 /*@
459    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
460    attempted by the nonlinear solver.
461 
462    Not Collective
463 
464    Input Parameter:
465 .  snes - SNES context
466 
467    Output Parameter:
468 .  nfails - number of unsuccessful steps attempted
469 
470    Notes:
471    This counter is reset to zero for each successive call to SNESSolve().
472 
473    Level: intermediate
474 
475 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
476 
477 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
478           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
479 @*/
480 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails)
481 {
482   PetscFunctionBegin;
483   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
484   PetscValidIntPointer(nfails,2);
485   *nfails = snes->numFailures;
486   PetscFunctionReturn(0);
487 }
488 
489 #undef __FUNCT__
490 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures"
491 /*@
492    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
493    attempted by the nonlinear solver before it gives up.
494 
495    Not Collective
496 
497    Input Parameters:
498 +  snes     - SNES context
499 -  maxFails - maximum of unsuccessful steps
500 
501    Level: intermediate
502 
503 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
504 
505 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
506           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
507 @*/
508 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
509 {
510   PetscFunctionBegin;
511   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
512   snes->maxFailures = maxFails;
513   PetscFunctionReturn(0);
514 }
515 
516 #undef __FUNCT__
517 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures"
518 /*@
519    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
520    attempted by the nonlinear solver before it gives up.
521 
522    Not Collective
523 
524    Input Parameter:
525 .  snes     - SNES context
526 
527    Output Parameter:
528 .  maxFails - maximum of unsuccessful steps
529 
530    Level: intermediate
531 
532 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
533 
534 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
535           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
536 
537 @*/
538 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
539 {
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
542   PetscValidIntPointer(maxFails,2);
543   *maxFails = snes->maxFailures;
544   PetscFunctionReturn(0);
545 }
546 
547 #undef __FUNCT__
548 #define __FUNCT__ "SNESGetNumberFunctionEvals"
549 /*@
550    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
551      done by SNES.
552 
553    Not Collective
554 
555    Input Parameter:
556 .  snes     - SNES context
557 
558    Output Parameter:
559 .  nfuncs - number of evaluations
560 
561    Level: intermediate
562 
563 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
564 
565 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures()
566 @*/
567 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
568 {
569   PetscFunctionBegin;
570   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
571   PetscValidIntPointer(nfuncs,2);
572   *nfuncs = snes->nfuncs;
573   PetscFunctionReturn(0);
574 }
575 
576 #undef __FUNCT__
577 #define __FUNCT__ "SNESGetLinearSolveFailures"
578 /*@
579    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
580    linear solvers.
581 
582    Not Collective
583 
584    Input Parameter:
585 .  snes - SNES context
586 
587    Output Parameter:
588 .  nfails - number of failed solves
589 
590    Notes:
591    This counter is reset to zero for each successive call to SNESSolve().
592 
593    Level: intermediate
594 
595 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
596 
597 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
598 @*/
599 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails)
600 {
601   PetscFunctionBegin;
602   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
603   PetscValidIntPointer(nfails,2);
604   *nfails = snes->numLinearSolveFailures;
605   PetscFunctionReturn(0);
606 }
607 
608 #undef __FUNCT__
609 #define __FUNCT__ "SNESSetMaxLinearSolveFailures"
610 /*@
611    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
612    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
613 
614    Collective on SNES
615 
616    Input Parameters:
617 +  snes     - SNES context
618 -  maxFails - maximum allowed linear solve failures
619 
620    Level: intermediate
621 
622    Notes: By default this is 0; that is SNES returns on the first failed linear solve
623 
624 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
625 
626 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
627 @*/
628 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
629 {
630   PetscFunctionBegin;
631   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
632   snes->maxLinearSolveFailures = maxFails;
633   PetscFunctionReturn(0);
634 }
635 
636 #undef __FUNCT__
637 #define __FUNCT__ "SNESGetMaxLinearSolveFailures"
638 /*@
639    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
640      are allowed before SNES terminates
641 
642    Not Collective
643 
644    Input Parameter:
645 .  snes     - SNES context
646 
647    Output Parameter:
648 .  maxFails - maximum of unsuccessful solves allowed
649 
650    Level: intermediate
651 
652    Notes: By default this is 1; that is SNES returns on the first failed linear solve
653 
654 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
655 
656 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
657 @*/
658 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
659 {
660   PetscFunctionBegin;
661   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
662   PetscValidIntPointer(maxFails,2);
663   *maxFails = snes->maxLinearSolveFailures;
664   PetscFunctionReturn(0);
665 }
666 
667 #undef __FUNCT__
668 #define __FUNCT__ "SNESGetLinearSolveIterations"
669 /*@
670    SNESGetLinearSolveIterations - Gets the total number of linear iterations
671    used by the nonlinear solver.
672 
673    Not Collective
674 
675    Input Parameter:
676 .  snes - SNES context
677 
678    Output Parameter:
679 .  lits - number of linear iterations
680 
681    Notes:
682    This counter is reset to zero for each successive call to SNESSolve().
683 
684    Level: intermediate
685 
686 .keywords: SNES, nonlinear, get, number, linear, iterations
687 
688 .seealso:  SNESGetIterationNumber(), SNESGetFunctionNorm()S, NESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures()
689 @*/
690 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveIterations(SNES snes,PetscInt* lits)
691 {
692   PetscFunctionBegin;
693   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
694   PetscValidIntPointer(lits,2);
695   *lits = snes->linear_its;
696   PetscFunctionReturn(0);
697 }
698 
699 #undef __FUNCT__
700 #define __FUNCT__ "SNESGetKSP"
701 /*@
702    SNESGetKSP - Returns the KSP context for a SNES solver.
703 
704    Not Collective, but if SNES object is parallel, then KSP object is parallel
705 
706    Input Parameter:
707 .  snes - the SNES context
708 
709    Output Parameter:
710 .  ksp - the KSP context
711 
712    Notes:
713    The user can then directly manipulate the KSP context to set various
714    options, etc.  Likewise, the user can then extract and manipulate the
715    PC contexts as well.
716 
717    Level: beginner
718 
719 .keywords: SNES, nonlinear, get, KSP, context
720 
721 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
722 @*/
723 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp)
724 {
725   PetscFunctionBegin;
726   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
727   PetscValidPointer(ksp,2);
728   *ksp = snes->ksp;
729   PetscFunctionReturn(0);
730 }
731 
732 #undef __FUNCT__
733 #define __FUNCT__ "SNESSetKSP"
734 /*@
735    SNESSetKSP - Sets a KSP context for the SNES object to use
736 
737    Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
738 
739    Input Parameters:
740 +  snes - the SNES context
741 -  ksp - the KSP context
742 
743    Notes:
744    The SNES object already has its KSP object, you can obtain with SNESGetKSP()
745    so this routine is rarely needed.
746 
747    The KSP object that is already in the SNES object has its reference count
748    decreased by one.
749 
750    Level: developer
751 
752 .keywords: SNES, nonlinear, get, KSP, context
753 
754 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
755 @*/
756 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetKSP(SNES snes,KSP ksp)
757 {
758   PetscErrorCode ierr;
759 
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
762   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
763   PetscCheckSameComm(snes,1,ksp,2);
764   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
765   if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);}
766   snes->ksp = ksp;
767   PetscFunctionReturn(0);
768 }
769 
770 #if 0
771 #undef __FUNCT__
772 #define __FUNCT__ "SNESPublish_Petsc"
773 static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
774 {
775   PetscFunctionBegin;
776   PetscFunctionReturn(0);
777 }
778 #endif
779 
780 /* -----------------------------------------------------------*/
781 #undef __FUNCT__
782 #define __FUNCT__ "SNESCreate"
783 /*@
784    SNESCreate - Creates a nonlinear solver context.
785 
786    Collective on MPI_Comm
787 
788    Input Parameters:
789 .  comm - MPI communicator
790 
791    Output Parameter:
792 .  outsnes - the new SNES context
793 
794    Options Database Keys:
795 +   -snes_mf - Activates default matrix-free Jacobian-vector products,
796                and no preconditioning matrix
797 .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
798                products, and a user-provided preconditioning matrix
799                as set by SNESSetJacobian()
800 -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
801 
802    Level: beginner
803 
804 .keywords: SNES, nonlinear, create, context
805 
806 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
807 
808 @*/
809 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes)
810 {
811   PetscErrorCode      ierr;
812   SNES                snes;
813   SNESKSPEW           *kctx;
814 
815   PetscFunctionBegin;
816   PetscValidPointer(outsnes,2);
817   *outsnes = PETSC_NULL;
818 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
819   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
820 #endif
821 
822   ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
823 
824   snes->ops->converged    = SNESDefaultConverged;
825   snes->max_its           = 50;
826   snes->max_funcs	  = 10000;
827   snes->norm		  = 0.0;
828   snes->rtol		  = 1.e-8;
829   snes->ttol              = 0.0;
830   snes->abstol		  = 1.e-50;
831   snes->xtol		  = 1.e-8;
832   snes->deltatol	  = 1.e-12;
833   snes->nfuncs            = 0;
834   snes->numFailures       = 0;
835   snes->maxFailures       = 1;
836   snes->linear_its        = 0;
837   snes->lagpreconditioner = 1;
838   snes->numbermonitors    = 0;
839   snes->data              = 0;
840   snes->setupcalled       = PETSC_FALSE;
841   snes->ksp_ewconv        = PETSC_FALSE;
842   snes->vwork             = 0;
843   snes->nwork             = 0;
844   snes->conv_hist_len     = 0;
845   snes->conv_hist_max     = 0;
846   snes->conv_hist         = PETSC_NULL;
847   snes->conv_hist_its     = PETSC_NULL;
848   snes->conv_hist_reset   = PETSC_TRUE;
849   snes->reason            = SNES_CONVERGED_ITERATING;
850 
851   snes->numLinearSolveFailures = 0;
852   snes->maxLinearSolveFailures = 1;
853 
854   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
855   ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr);
856   snes->kspconvctx  = (void*)kctx;
857   kctx->version     = 2;
858   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
859                              this was too large for some test cases */
860   kctx->rtol_last   = 0;
861   kctx->rtol_max    = .9;
862   kctx->gamma       = 1.0;
863   kctx->alpha       = .5*(1.0 + sqrt(5.0));
864   kctx->alpha2      = kctx->alpha;
865   kctx->threshold   = .1;
866   kctx->lresid_last = 0;
867   kctx->norm_last   = 0;
868 
869   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
870   ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
871 
872   *outsnes = snes;
873   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
874   PetscFunctionReturn(0);
875 }
876 
877 #undef __FUNCT__
878 #define __FUNCT__ "SNESSetFunction"
879 /*@C
880    SNESSetFunction - Sets the function evaluation routine and function
881    vector for use by the SNES routines in solving systems of nonlinear
882    equations.
883 
884    Collective on SNES
885 
886    Input Parameters:
887 +  snes - the SNES context
888 .  r - vector to store function value
889 .  func - function evaluation routine
890 -  ctx - [optional] user-defined context for private data for the
891          function evaluation routine (may be PETSC_NULL)
892 
893    Calling sequence of func:
894 $    func (SNES snes,Vec x,Vec f,void *ctx);
895 
896 .  f - function vector
897 -  ctx - optional user-defined function context
898 
899    Notes:
900    The Newton-like methods typically solve linear systems of the form
901 $      f'(x) x = -f(x),
902    where f'(x) denotes the Jacobian matrix and f(x) is the function.
903 
904    Level: beginner
905 
906 .keywords: SNES, nonlinear, set, function
907 
908 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
909 @*/
910 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
911 {
912   PetscErrorCode ierr;
913   PetscFunctionBegin;
914   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
915   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
916   PetscCheckSameComm(snes,1,r,2);
917   ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr);
918   if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); }
919   snes->ops->computefunction = func;
920   snes->vec_func             = r;
921   snes->funP                 = ctx;
922   PetscFunctionReturn(0);
923 }
924 
925 /* --------------------------------------------------------------- */
926 #undef __FUNCT__
927 #define __FUNCT__ "SNESGetRhs"
928 /*@C
929    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
930    it assumes a zero right hand side.
931 
932    Collective on SNES
933 
934    Input Parameter:
935 .  snes - the SNES context
936 
937    Output Parameter:
938 .  rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null
939 
940    Level: intermediate
941 
942 .keywords: SNES, nonlinear, get, function, right hand side
943 
944 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
945 @*/
946 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs)
947 {
948   PetscFunctionBegin;
949   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
950   PetscValidPointer(rhs,2);
951   *rhs = snes->vec_rhs;
952   PetscFunctionReturn(0);
953 }
954 
955 #undef __FUNCT__
956 #define __FUNCT__ "SNESComputeFunction"
957 /*@
958    SNESComputeFunction - Calls the function that has been set with
959                          SNESSetFunction().
960 
961    Collective on SNES
962 
963    Input Parameters:
964 +  snes - the SNES context
965 -  x - input vector
966 
967    Output Parameter:
968 .  y - function vector, as set by SNESSetFunction()
969 
970    Notes:
971    SNESComputeFunction() is typically used within nonlinear solvers
972    implementations, so most users would not generally call this routine
973    themselves.
974 
975    Level: developer
976 
977 .keywords: SNES, nonlinear, compute, function
978 
979 .seealso: SNESSetFunction(), SNESGetFunction()
980 @*/
981 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y)
982 {
983   PetscErrorCode ierr;
984 
985   PetscFunctionBegin;
986   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
987   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
988   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
989   PetscCheckSameComm(snes,1,x,2);
990   PetscCheckSameComm(snes,1,y,3);
991 
992   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
993   if (snes->ops->computefunction) {
994     PetscStackPush("SNES user function");
995     CHKMEMQ;
996     ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);
997     CHKMEMQ;
998     PetscStackPop;
999     if (PetscExceptionValue(ierr)) {
1000       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
1001     }
1002     CHKERRQ(ierr);
1003   } else if (snes->vec_rhs) {
1004     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
1005   } else {
1006     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
1007   }
1008   if (snes->vec_rhs) {
1009     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
1010   }
1011   snes->nfuncs++;
1012   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 #undef __FUNCT__
1017 #define __FUNCT__ "SNESComputeJacobian"
1018 /*@
1019    SNESComputeJacobian - Computes the Jacobian matrix that has been
1020    set with SNESSetJacobian().
1021 
1022    Collective on SNES and Mat
1023 
1024    Input Parameters:
1025 +  snes - the SNES context
1026 -  x - input vector
1027 
1028    Output Parameters:
1029 +  A - Jacobian matrix
1030 .  B - optional preconditioning matrix
1031 -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1032 
1033    Notes:
1034    Most users should not need to explicitly call this routine, as it
1035    is used internally within the nonlinear solvers.
1036 
1037    See KSPSetOperators() for important information about setting the
1038    flag parameter.
1039 
1040    Level: developer
1041 
1042 .keywords: SNES, compute, Jacobian, matrix
1043 
1044 .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure
1045 @*/
1046 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1047 {
1048   PetscErrorCode ierr;
1049 
1050   PetscFunctionBegin;
1051   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1052   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
1053   PetscValidPointer(flg,5);
1054   PetscCheckSameComm(snes,1,X,2);
1055   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1056   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1057   *flg = DIFFERENT_NONZERO_PATTERN;
1058   PetscStackPush("SNES user Jacobian function");
1059   CHKMEMQ;
1060   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1061   CHKMEMQ;
1062   PetscStackPop;
1063   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1064 
1065   if (snes->lagpreconditioner == -1) {
1066     *flg = SAME_PRECONDITIONER;
1067     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1068   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1069     *flg = SAME_PRECONDITIONER;
1070     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1071   }
1072 
1073   /* make sure user returned a correct Jacobian and preconditioner */
1074   /* PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
1075     PetscValidHeaderSpecific(*B,MAT_COOKIE,4);   */
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 #undef __FUNCT__
1080 #define __FUNCT__ "SNESSetJacobian"
1081 /*@C
1082    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1083    location to store the matrix.
1084 
1085    Collective on SNES and Mat
1086 
1087    Input Parameters:
1088 +  snes - the SNES context
1089 .  A - Jacobian matrix
1090 .  B - preconditioner matrix (usually same as the Jacobian)
1091 .  func - Jacobian evaluation routine
1092 -  ctx - [optional] user-defined context for private data for the
1093          Jacobian evaluation routine (may be PETSC_NULL)
1094 
1095    Calling sequence of func:
1096 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1097 
1098 +  x - input vector
1099 .  A - Jacobian matrix
1100 .  B - preconditioner matrix, usually the same as A
1101 .  flag - flag indicating information about the preconditioner matrix
1102    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1103 -  ctx - [optional] user-defined Jacobian context
1104 
1105    Notes:
1106    See KSPSetOperators() for important information about setting the flag
1107    output parameter in the routine func().  Be sure to read this information!
1108 
1109    The routine func() takes Mat * as the matrix arguments rather than Mat.
1110    This allows the Jacobian evaluation routine to replace A and/or B with a
1111    completely new new matrix structure (not just different matrix elements)
1112    when appropriate, for instance, if the nonzero structure is changing
1113    throughout the global iterations.
1114 
1115    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
1116    each matrix.
1117 
1118    Level: beginner
1119 
1120 .keywords: SNES, nonlinear, set, Jacobian, matrix
1121 
1122 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
1123 @*/
1124 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1125 {
1126   PetscErrorCode ierr;
1127 
1128   PetscFunctionBegin;
1129   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1130   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
1131   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
1132   if (A) PetscCheckSameComm(snes,1,A,2);
1133   if (B) PetscCheckSameComm(snes,1,B,2);
1134   if (func) snes->ops->computejacobian = func;
1135   if (ctx)  snes->jacP                 = ctx;
1136   if (A) {
1137     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1138     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1139     snes->jacobian = A;
1140   }
1141   if (B) {
1142     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1143     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1144     snes->jacobian_pre = B;
1145   }
1146   PetscFunctionReturn(0);
1147 }
1148 
1149 #undef __FUNCT__
1150 #define __FUNCT__ "SNESGetJacobian"
1151 /*@C
1152    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1153    provided context for evaluating the Jacobian.
1154 
1155    Not Collective, but Mat object will be parallel if SNES object is
1156 
1157    Input Parameter:
1158 .  snes - the nonlinear solver context
1159 
1160    Output Parameters:
1161 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1162 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1163 .  func - location to put Jacobian function (or PETSC_NULL)
1164 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1165 
1166    Level: advanced
1167 
1168 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1169 @*/
1170 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1171 {
1172   PetscFunctionBegin;
1173   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1174   if (A)    *A    = snes->jacobian;
1175   if (B)    *B    = snes->jacobian_pre;
1176   if (func) *func = snes->ops->computejacobian;
1177   if (ctx)  *ctx  = snes->jacP;
1178   PetscFunctionReturn(0);
1179 }
1180 
1181 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1182 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
1183 
1184 #undef __FUNCT__
1185 #define __FUNCT__ "SNESSetUp"
1186 /*@
1187    SNESSetUp - Sets up the internal data structures for the later use
1188    of a nonlinear solver.
1189 
1190    Collective on SNES
1191 
1192    Input Parameters:
1193 .  snes - the SNES context
1194 
1195    Notes:
1196    For basic use of the SNES solvers the user need not explicitly call
1197    SNESSetUp(), since these actions will automatically occur during
1198    the call to SNESSolve().  However, if one wishes to control this
1199    phase separately, SNESSetUp() should be called after SNESCreate()
1200    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1201 
1202    Level: advanced
1203 
1204 .keywords: SNES, nonlinear, setup
1205 
1206 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1207 @*/
1208 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
1209 {
1210   PetscErrorCode ierr;
1211   PetscTruth     flg;
1212 
1213   PetscFunctionBegin;
1214   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1215   if (snes->setupcalled) PetscFunctionReturn(0);
1216 
1217   if (!((PetscObject)snes)->type_name) {
1218     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
1219   }
1220 
1221   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1222   /*
1223       This version replaces the user provided Jacobian matrix with a
1224       matrix-free version but still employs the user-provided preconditioner matrix
1225   */
1226   if (flg) {
1227     Mat J;
1228     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1229     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1230     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1231     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1232     ierr = MatDestroy(J);CHKERRQ(ierr);
1233   }
1234 
1235 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT)
1236   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
1237   if (flg) {
1238     Mat J;
1239     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1240     ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr);
1241     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
1242     ierr = MatDestroy(J);CHKERRQ(ierr);
1243   }
1244 #endif
1245 
1246   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1247   /*
1248       This version replaces both the user-provided Jacobian and the user-
1249       provided preconditioner matrix with the default matrix free version.
1250    */
1251   if (flg) {
1252     Mat  J;
1253     KSP ksp;
1254     PC   pc;
1255     /* create and set matrix-free operator */
1256     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1257     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1258     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
1259     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
1260     ierr = MatDestroy(J);CHKERRQ(ierr);
1261     /* force no preconditioner */
1262     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1263     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1264     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1265     if (!flg) {
1266       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
1267       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1268     }
1269   }
1270 
1271   if (!snes->vec_func && !snes->vec_rhs) {
1272     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1273   }
1274   if (!snes->ops->computefunction && !snes->vec_rhs) {
1275     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1276   }
1277   if (!snes->jacobian) {
1278     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1279   }
1280   if (snes->vec_func == snes->vec_sol) {
1281     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1282   }
1283 
1284   if (snes->ops->setup) {
1285     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1286   }
1287   snes->setupcalled = PETSC_TRUE;
1288   PetscFunctionReturn(0);
1289 }
1290 
1291 #undef __FUNCT__
1292 #define __FUNCT__ "SNESDestroy"
1293 /*@
1294    SNESDestroy - Destroys the nonlinear solver context that was created
1295    with SNESCreate().
1296 
1297    Collective on SNES
1298 
1299    Input Parameter:
1300 .  snes - the SNES context
1301 
1302    Level: beginner
1303 
1304 .keywords: SNES, nonlinear, destroy
1305 
1306 .seealso: SNESCreate(), SNESSolve()
1307 @*/
1308 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
1309 {
1310   PetscErrorCode ierr;
1311 
1312   PetscFunctionBegin;
1313   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1314   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1315 
1316   /* if memory was published with AMS then destroy it */
1317   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1318   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
1319 
1320   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
1321   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
1322   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
1323   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
1324   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
1325   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1326   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1327   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1328   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
1329   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1330   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
1331   PetscFunctionReturn(0);
1332 }
1333 
1334 /* ----------- Routines to set solver parameters ---------- */
1335 
1336 #undef __FUNCT__
1337 #define __FUNCT__ "SNESSetLagPreconditioner"
1338 /*@
1339    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1340 
1341    Collective on SNES
1342 
1343    Input Parameters:
1344 +  snes - the SNES context
1345 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1346          the Jacobian is built etc.
1347 
1348    Options Database Keys:
1349 .    -snes_lag_preconditioner <lag>
1350 
1351    Notes:
1352    The default is 1
1353    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1354    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1355 
1356    Level: intermediate
1357 
1358 .keywords: SNES, nonlinear, set, convergence, tolerances
1359 
1360 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner()
1361 
1362 @*/
1363 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1364 {
1365   PetscFunctionBegin;
1366   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1367   if (lag < -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -1, 1 or greater");
1368   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1369   snes->lagpreconditioner = lag;
1370   PetscFunctionReturn(0);
1371 }
1372 
1373 #undef __FUNCT__
1374 #define __FUNCT__ "SNESGetLagPreconditioner"
1375 /*@
1376    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1377 
1378    Collective on SNES
1379 
1380    Input Parameter:
1381 .  snes - the SNES context
1382 
1383    Output Parameter:
1384 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1385          the Jacobian is built etc.
1386 
1387    Options Database Keys:
1388 .    -snes_lag_preconditioner <lag>
1389 
1390    Notes:
1391    The default is 1
1392    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1393 
1394    Level: intermediate
1395 
1396 .keywords: SNES, nonlinear, set, convergence, tolerances
1397 
1398 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1399 
1400 @*/
1401 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1402 {
1403   PetscFunctionBegin;
1404   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1405   *lag = snes->lagpreconditioner;
1406   PetscFunctionReturn(0);
1407 }
1408 
1409 #undef __FUNCT__
1410 #define __FUNCT__ "SNESSetTolerances"
1411 /*@
1412    SNESSetTolerances - Sets various parameters used in convergence tests.
1413 
1414    Collective on SNES
1415 
1416    Input Parameters:
1417 +  snes - the SNES context
1418 .  abstol - absolute convergence tolerance
1419 .  rtol - relative convergence tolerance
1420 .  stol -  convergence tolerance in terms of the norm
1421            of the change in the solution between steps
1422 .  maxit - maximum number of iterations
1423 -  maxf - maximum number of function evaluations
1424 
1425    Options Database Keys:
1426 +    -snes_atol <abstol> - Sets abstol
1427 .    -snes_rtol <rtol> - Sets rtol
1428 .    -snes_stol <stol> - Sets stol
1429 .    -snes_max_it <maxit> - Sets maxit
1430 -    -snes_max_funcs <maxf> - Sets maxf
1431 
1432    Notes:
1433    The default maximum number of iterations is 50.
1434    The default maximum number of function evaluations is 1000.
1435 
1436    Level: intermediate
1437 
1438 .keywords: SNES, nonlinear, set, convergence, tolerances
1439 
1440 .seealso: SNESSetTrustRegionTolerance()
1441 @*/
1442 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
1443 {
1444   PetscFunctionBegin;
1445   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1446   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1447   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1448   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1449   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1450   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
1451   PetscFunctionReturn(0);
1452 }
1453 
1454 #undef __FUNCT__
1455 #define __FUNCT__ "SNESGetTolerances"
1456 /*@
1457    SNESGetTolerances - Gets various parameters used in convergence tests.
1458 
1459    Not Collective
1460 
1461    Input Parameters:
1462 +  snes - the SNES context
1463 .  atol - absolute convergence tolerance
1464 .  rtol - relative convergence tolerance
1465 .  stol -  convergence tolerance in terms of the norm
1466            of the change in the solution between steps
1467 .  maxit - maximum number of iterations
1468 -  maxf - maximum number of function evaluations
1469 
1470    Notes:
1471    The user can specify PETSC_NULL for any parameter that is not needed.
1472 
1473    Level: intermediate
1474 
1475 .keywords: SNES, nonlinear, get, convergence, tolerances
1476 
1477 .seealso: SNESSetTolerances()
1478 @*/
1479 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
1480 {
1481   PetscFunctionBegin;
1482   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1483   if (atol)  *atol  = snes->abstol;
1484   if (rtol)  *rtol  = snes->rtol;
1485   if (stol)  *stol  = snes->xtol;
1486   if (maxit) *maxit = snes->max_its;
1487   if (maxf)  *maxf  = snes->max_funcs;
1488   PetscFunctionReturn(0);
1489 }
1490 
1491 #undef __FUNCT__
1492 #define __FUNCT__ "SNESSetTrustRegionTolerance"
1493 /*@
1494    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
1495 
1496    Collective on SNES
1497 
1498    Input Parameters:
1499 +  snes - the SNES context
1500 -  tol - tolerance
1501 
1502    Options Database Key:
1503 .  -snes_trtol <tol> - Sets tol
1504 
1505    Level: intermediate
1506 
1507 .keywords: SNES, nonlinear, set, trust region, tolerance
1508 
1509 .seealso: SNESSetTolerances()
1510 @*/
1511 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
1512 {
1513   PetscFunctionBegin;
1514   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1515   snes->deltatol = tol;
1516   PetscFunctionReturn(0);
1517 }
1518 
1519 /*
1520    Duplicate the lg monitors for SNES from KSP; for some reason with
1521    dynamic libraries things don't work under Sun4 if we just use
1522    macros instead of functions
1523 */
1524 #undef __FUNCT__
1525 #define __FUNCT__ "SNESMonitorLG"
1526 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1527 {
1528   PetscErrorCode ierr;
1529 
1530   PetscFunctionBegin;
1531   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1532   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1533   PetscFunctionReturn(0);
1534 }
1535 
1536 #undef __FUNCT__
1537 #define __FUNCT__ "SNESMonitorLGCreate"
1538 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1539 {
1540   PetscErrorCode ierr;
1541 
1542   PetscFunctionBegin;
1543   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1544   PetscFunctionReturn(0);
1545 }
1546 
1547 #undef __FUNCT__
1548 #define __FUNCT__ "SNESMonitorLGDestroy"
1549 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1550 {
1551   PetscErrorCode ierr;
1552 
1553   PetscFunctionBegin;
1554   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1555   PetscFunctionReturn(0);
1556 }
1557 
1558 /* ------------ Routines to set performance monitoring options ----------- */
1559 
1560 #undef __FUNCT__
1561 #define __FUNCT__ "SNESMonitorSet"
1562 /*@C
1563    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
1564    iteration of the nonlinear solver to display the iteration's
1565    progress.
1566 
1567    Collective on SNES
1568 
1569    Input Parameters:
1570 +  snes - the SNES context
1571 .  func - monitoring routine
1572 .  mctx - [optional] user-defined context for private data for the
1573           monitor routine (use PETSC_NULL if no context is desired)
1574 -  monitordestroy - [optional] routine that frees monitor context
1575           (may be PETSC_NULL)
1576 
1577    Calling sequence of func:
1578 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1579 
1580 +    snes - the SNES context
1581 .    its - iteration number
1582 .    norm - 2-norm function value (may be estimated)
1583 -    mctx - [optional] monitoring context
1584 
1585    Options Database Keys:
1586 +    -snes_monitor        - sets SNESMonitorDefault()
1587 .    -snes_monitor_draw    - sets line graph monitor,
1588                             uses SNESMonitorLGCreate()
1589 _    -snes_monitor_cancel - cancels all monitors that have
1590                             been hardwired into a code by
1591                             calls to SNESMonitorSet(), but
1592                             does not cancel those set via
1593                             the options database.
1594 
1595    Notes:
1596    Several different monitoring routines may be set by calling
1597    SNESMonitorSet() multiple times; all will be called in the
1598    order in which they were set.
1599 
1600    Fortran notes: Only a single monitor function can be set for each SNES object
1601 
1602    Level: intermediate
1603 
1604 .keywords: SNES, nonlinear, set, monitor
1605 
1606 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
1607 @*/
1608 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
1609 {
1610   PetscInt i;
1611 
1612   PetscFunctionBegin;
1613   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1614   if (snes->numbermonitors >= MAXSNESMONITORS) {
1615     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1616   }
1617   for (i=0; i<snes->numbermonitors;i++) {
1618     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1619 
1620     /* check if both default monitors that share common ASCII viewer */
1621     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1622       if (mctx && snes->monitorcontext[i]) {
1623         PetscErrorCode          ierr;
1624         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1625         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1626         if (viewer1->viewer == viewer2->viewer) {
1627           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1628           PetscFunctionReturn(0);
1629         }
1630       }
1631     }
1632   }
1633   snes->monitor[snes->numbermonitors]           = monitor;
1634   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1635   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
1636   PetscFunctionReturn(0);
1637 }
1638 
1639 #undef __FUNCT__
1640 #define __FUNCT__ "SNESMonitorCancel"
1641 /*@C
1642    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
1643 
1644    Collective on SNES
1645 
1646    Input Parameters:
1647 .  snes - the SNES context
1648 
1649    Options Database Key:
1650 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1651     into a code by calls to SNESMonitorSet(), but does not cancel those
1652     set via the options database
1653 
1654    Notes:
1655    There is no way to clear one specific monitor from a SNES object.
1656 
1657    Level: intermediate
1658 
1659 .keywords: SNES, nonlinear, set, monitor
1660 
1661 .seealso: SNESMonitorDefault(), SNESMonitorSet()
1662 @*/
1663 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
1664 {
1665   PetscErrorCode ierr;
1666   PetscInt       i;
1667 
1668   PetscFunctionBegin;
1669   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1670   for (i=0; i<snes->numbermonitors; i++) {
1671     if (snes->monitordestroy[i]) {
1672       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1673     }
1674   }
1675   snes->numbermonitors = 0;
1676   PetscFunctionReturn(0);
1677 }
1678 
1679 #undef __FUNCT__
1680 #define __FUNCT__ "SNESSetConvergenceTest"
1681 /*@C
1682    SNESSetConvergenceTest - Sets the function that is to be used
1683    to test for convergence of the nonlinear iterative solution.
1684 
1685    Collective on SNES
1686 
1687    Input Parameters:
1688 +  snes - the SNES context
1689 .  func - routine to test for convergence
1690 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
1691 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
1692 
1693    Calling sequence of func:
1694 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1695 
1696 +    snes - the SNES context
1697 .    it - current iteration (0 is the first and is before any Newton step)
1698 .    cctx - [optional] convergence context
1699 .    reason - reason for convergence/divergence
1700 .    xnorm - 2-norm of current iterate
1701 .    gnorm - 2-norm of current step
1702 -    f - 2-norm of function
1703 
1704    Level: advanced
1705 
1706 .keywords: SNES, nonlinear, set, convergence, test
1707 
1708 .seealso: SNESDefaultConverged(), SNESSkipConverged()
1709 @*/
1710 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1711 {
1712   PetscErrorCode ierr;
1713 
1714   PetscFunctionBegin;
1715   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1716   if (!func) func = SNESSkipConverged;
1717   if (snes->ops->convergeddestroy) {
1718     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
1719   }
1720   snes->ops->converged        = func;
1721   snes->ops->convergeddestroy = destroy;
1722   snes->cnvP                  = cctx;
1723   PetscFunctionReturn(0);
1724 }
1725 
1726 #undef __FUNCT__
1727 #define __FUNCT__ "SNESGetConvergedReason"
1728 /*@
1729    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1730 
1731    Not Collective
1732 
1733    Input Parameter:
1734 .  snes - the SNES context
1735 
1736    Output Parameter:
1737 .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1738             manual pages for the individual convergence tests for complete lists
1739 
1740    Level: intermediate
1741 
1742    Notes: Can only be called after the call the SNESSolve() is complete.
1743 
1744 .keywords: SNES, nonlinear, set, convergence, test
1745 
1746 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
1747 @*/
1748 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1749 {
1750   PetscFunctionBegin;
1751   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1752   PetscValidPointer(reason,2);
1753   *reason = snes->reason;
1754   PetscFunctionReturn(0);
1755 }
1756 
1757 #undef __FUNCT__
1758 #define __FUNCT__ "SNESSetConvergenceHistory"
1759 /*@
1760    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1761 
1762    Collective on SNES
1763 
1764    Input Parameters:
1765 +  snes - iterative context obtained from SNESCreate()
1766 .  a   - array to hold history
1767 .  its - integer array holds the number of linear iterations for each solve.
1768 .  na  - size of a and its
1769 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1770            else it continues storing new values for new nonlinear solves after the old ones
1771 
1772    Notes:
1773    If set, this array will contain the function norms computed
1774    at each step.
1775 
1776    This routine is useful, e.g., when running a code for purposes
1777    of accurate performance monitoring, when no I/O should be done
1778    during the section of code that is being timed.
1779 
1780    Level: intermediate
1781 
1782 .keywords: SNES, set, convergence, history
1783 
1784 .seealso: SNESGetConvergenceHistory()
1785 
1786 @*/
1787 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset)
1788 {
1789   PetscFunctionBegin;
1790   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1791   if (na)  PetscValidScalarPointer(a,2);
1792   if (its) PetscValidIntPointer(its,3);
1793   snes->conv_hist       = a;
1794   snes->conv_hist_its   = its;
1795   snes->conv_hist_max   = na;
1796   snes->conv_hist_len   = 0;
1797   snes->conv_hist_reset = reset;
1798   PetscFunctionReturn(0);
1799 }
1800 
1801 #undef __FUNCT__
1802 #define __FUNCT__ "SNESGetConvergenceHistory"
1803 /*@C
1804    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1805 
1806    Collective on SNES
1807 
1808    Input Parameter:
1809 .  snes - iterative context obtained from SNESCreate()
1810 
1811    Output Parameters:
1812 .  a   - array to hold history
1813 .  its - integer array holds the number of linear iterations (or
1814          negative if not converged) for each solve.
1815 -  na  - size of a and its
1816 
1817    Notes:
1818     The calling sequence for this routine in Fortran is
1819 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1820 
1821    This routine is useful, e.g., when running a code for purposes
1822    of accurate performance monitoring, when no I/O should be done
1823    during the section of code that is being timed.
1824 
1825    Level: intermediate
1826 
1827 .keywords: SNES, get, convergence, history
1828 
1829 .seealso: SNESSetConvergencHistory()
1830 
1831 @*/
1832 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1833 {
1834   PetscFunctionBegin;
1835   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1836   if (a)   *a   = snes->conv_hist;
1837   if (its) *its = snes->conv_hist_its;
1838   if (na)  *na  = snes->conv_hist_len;
1839   PetscFunctionReturn(0);
1840 }
1841 
1842 #undef __FUNCT__
1843 #define __FUNCT__ "SNESSetUpdate"
1844 /*@C
1845   SNESSetUpdate - Sets the general-purpose update function called
1846   at the beginning o every iteration of the nonlinear solve. Specifically
1847   it is called just before the Jacobian is "evaluated".
1848 
1849   Collective on SNES
1850 
1851   Input Parameters:
1852 . snes - The nonlinear solver context
1853 . func - The function
1854 
1855   Calling sequence of func:
1856 . func (SNES snes, PetscInt step);
1857 
1858 . step - The current step of the iteration
1859 
1860   Level: intermediate
1861 
1862 .keywords: SNES, update
1863 
1864 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
1865 @*/
1866 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
1867 {
1868   PetscFunctionBegin;
1869   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
1870   snes->ops->update = func;
1871   PetscFunctionReturn(0);
1872 }
1873 
1874 #undef __FUNCT__
1875 #define __FUNCT__ "SNESDefaultUpdate"
1876 /*@
1877   SNESDefaultUpdate - The default update function which does nothing.
1878 
1879   Not collective
1880 
1881   Input Parameters:
1882 . snes - The nonlinear solver context
1883 . step - The current step of the iteration
1884 
1885   Level: intermediate
1886 
1887 .keywords: SNES, update
1888 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
1889 @*/
1890 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
1891 {
1892   PetscFunctionBegin;
1893   PetscFunctionReturn(0);
1894 }
1895 
1896 #undef __FUNCT__
1897 #define __FUNCT__ "SNESScaleStep_Private"
1898 /*
1899    SNESScaleStep_Private - Scales a step so that its length is less than the
1900    positive parameter delta.
1901 
1902     Input Parameters:
1903 +   snes - the SNES context
1904 .   y - approximate solution of linear system
1905 .   fnorm - 2-norm of current function
1906 -   delta - trust region size
1907 
1908     Output Parameters:
1909 +   gpnorm - predicted function norm at the new point, assuming local
1910     linearization.  The value is zero if the step lies within the trust
1911     region, and exceeds zero otherwise.
1912 -   ynorm - 2-norm of the step
1913 
1914     Note:
1915     For non-trust region methods such as SNESLS, the parameter delta
1916     is set to be the maximum allowable step size.
1917 
1918 .keywords: SNES, nonlinear, scale, step
1919 */
1920 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
1921 {
1922   PetscReal      nrm;
1923   PetscScalar    cnorm;
1924   PetscErrorCode ierr;
1925 
1926   PetscFunctionBegin;
1927   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1928   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1929   PetscCheckSameComm(snes,1,y,2);
1930 
1931   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1932   if (nrm > *delta) {
1933      nrm = *delta/nrm;
1934      *gpnorm = (1.0 - nrm)*(*fnorm);
1935      cnorm = nrm;
1936      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
1937      *ynorm = *delta;
1938   } else {
1939      *gpnorm = 0.0;
1940      *ynorm = nrm;
1941   }
1942   PetscFunctionReturn(0);
1943 }
1944 
1945 #undef __FUNCT__
1946 #define __FUNCT__ "SNESSolve"
1947 /*@C
1948    SNESSolve - Solves a nonlinear system F(x) = b.
1949    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
1950 
1951    Collective on SNES
1952 
1953    Input Parameters:
1954 +  snes - the SNES context
1955 .  b - the constant part of the equation, or PETSC_NULL to use zero.
1956 -  x - the solution vector.
1957 
1958    Notes:
1959    The user should initialize the vector,x, with the initial guess
1960    for the nonlinear solve prior to calling SNESSolve.  In particular,
1961    to employ an initial guess of zero, the user should explicitly set
1962    this vector to zero by calling VecSet().
1963 
1964    Level: beginner
1965 
1966 .keywords: SNES, nonlinear, solve
1967 
1968 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
1969 @*/
1970 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
1971 {
1972   PetscErrorCode ierr;
1973   PetscTruth     flg;
1974   char           filename[PETSC_MAX_PATH_LEN];
1975   PetscViewer    viewer;
1976 
1977   PetscFunctionBegin;
1978   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1979   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
1980   PetscCheckSameComm(snes,1,x,3);
1981   if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2);
1982   if (b) PetscCheckSameComm(snes,1,b,2);
1983 
1984   /* set solution vector */
1985   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
1986   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
1987   snes->vec_sol = x;
1988   /* set afine vector if provided */
1989   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
1990   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
1991   snes->vec_rhs = b;
1992 
1993   if (!snes->vec_func && snes->vec_rhs) {
1994     ierr = VecDuplicate(b, &snes->vec_func); CHKERRQ(ierr);
1995   }
1996 
1997   ierr = SNESSetUp(snes);CHKERRQ(ierr);
1998 
1999   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
2000   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2001 
2002   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2003   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
2004   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2005   if (snes->domainerror){
2006     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
2007     snes->domainerror = PETSC_FALSE;
2008   }
2009 
2010   if (!snes->reason) {
2011     SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
2012   }
2013 
2014   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2015   if (flg && !PetscPreLoadingOn) {
2016     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2017     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2018     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2019   }
2020 
2021   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
2022   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
2023   if (snes->printreason) {
2024     if (snes->reason > 0) {
2025       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2026     } else {
2027       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
2028     }
2029   }
2030 
2031   PetscFunctionReturn(0);
2032 }
2033 
2034 /* --------- Internal routines for SNES Package --------- */
2035 
2036 #undef __FUNCT__
2037 #define __FUNCT__ "SNESSetType"
2038 /*@C
2039    SNESSetType - Sets the method for the nonlinear solver.
2040 
2041    Collective on SNES
2042 
2043    Input Parameters:
2044 +  snes - the SNES context
2045 -  type - a known method
2046 
2047    Options Database Key:
2048 .  -snes_type <type> - Sets the method; use -help for a list
2049    of available methods (for instance, ls or tr)
2050 
2051    Notes:
2052    See "petsc/include/petscsnes.h" for available methods (for instance)
2053 +    SNESLS - Newton's method with line search
2054      (systems of nonlinear equations)
2055 .    SNESTR - Newton's method with trust region
2056      (systems of nonlinear equations)
2057 
2058   Normally, it is best to use the SNESSetFromOptions() command and then
2059   set the SNES solver type from the options database rather than by using
2060   this routine.  Using the options database provides the user with
2061   maximum flexibility in evaluating the many nonlinear solvers.
2062   The SNESSetType() routine is provided for those situations where it
2063   is necessary to set the nonlinear solver independently of the command
2064   line or options database.  This might be the case, for example, when
2065   the choice of solver changes during the execution of the program,
2066   and the user's application is taking responsibility for choosing the
2067   appropriate method.
2068 
2069   Level: intermediate
2070 
2071 .keywords: SNES, set, type
2072 
2073 .seealso: SNESType, SNESCreate()
2074 
2075 @*/
2076 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type)
2077 {
2078   PetscErrorCode ierr,(*r)(SNES);
2079   PetscTruth     match;
2080 
2081   PetscFunctionBegin;
2082   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2083   PetscValidCharPointer(type,2);
2084 
2085   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
2086   if (match) PetscFunctionReturn(0);
2087 
2088   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2089   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
2090   /* Destroy the previous private SNES context */
2091   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
2092   /* Reinitialize function pointers in SNESOps structure */
2093   snes->ops->setup          = 0;
2094   snes->ops->solve          = 0;
2095   snes->ops->view           = 0;
2096   snes->ops->setfromoptions = 0;
2097   snes->ops->destroy        = 0;
2098   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
2099   snes->setupcalled = PETSC_FALSE;
2100   ierr = (*r)(snes);CHKERRQ(ierr);
2101   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
2102   PetscFunctionReturn(0);
2103 }
2104 
2105 
2106 /* --------------------------------------------------------------------- */
2107 #undef __FUNCT__
2108 #define __FUNCT__ "SNESRegisterDestroy"
2109 /*@
2110    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2111    registered by SNESRegisterDynamic().
2112 
2113    Not Collective
2114 
2115    Level: advanced
2116 
2117 .keywords: SNES, nonlinear, register, destroy
2118 
2119 .seealso: SNESRegisterAll(), SNESRegisterAll()
2120 @*/
2121 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
2122 {
2123   PetscErrorCode ierr;
2124 
2125   PetscFunctionBegin;
2126   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
2127   SNESRegisterAllCalled = PETSC_FALSE;
2128   PetscFunctionReturn(0);
2129 }
2130 
2131 #undef __FUNCT__
2132 #define __FUNCT__ "SNESGetType"
2133 /*@C
2134    SNESGetType - Gets the SNES method type and name (as a string).
2135 
2136    Not Collective
2137 
2138    Input Parameter:
2139 .  snes - nonlinear solver context
2140 
2141    Output Parameter:
2142 .  type - SNES method (a character string)
2143 
2144    Level: intermediate
2145 
2146 .keywords: SNES, nonlinear, get, type, name
2147 @*/
2148 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type)
2149 {
2150   PetscFunctionBegin;
2151   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2152   PetscValidPointer(type,2);
2153   *type = ((PetscObject)snes)->type_name;
2154   PetscFunctionReturn(0);
2155 }
2156 
2157 #undef __FUNCT__
2158 #define __FUNCT__ "SNESGetSolution"
2159 /*@
2160    SNESGetSolution - Returns the vector where the approximate solution is
2161    stored.
2162 
2163    Not Collective, but Vec is parallel if SNES is parallel
2164 
2165    Input Parameter:
2166 .  snes - the SNES context
2167 
2168    Output Parameter:
2169 .  x - the solution
2170 
2171    Level: intermediate
2172 
2173 .keywords: SNES, nonlinear, get, solution
2174 
2175 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
2176 @*/
2177 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
2178 {
2179   PetscFunctionBegin;
2180   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2181   PetscValidPointer(x,2);
2182   *x = snes->vec_sol;
2183   PetscFunctionReturn(0);
2184 }
2185 
2186 #undef __FUNCT__
2187 #define __FUNCT__ "SNESGetSolutionUpdate"
2188 /*@
2189    SNESGetSolutionUpdate - Returns the vector where the solution update is
2190    stored.
2191 
2192    Not Collective, but Vec is parallel if SNES is parallel
2193 
2194    Input Parameter:
2195 .  snes - the SNES context
2196 
2197    Output Parameter:
2198 .  x - the solution update
2199 
2200    Level: advanced
2201 
2202 .keywords: SNES, nonlinear, get, solution, update
2203 
2204 .seealso: SNESGetSolution(), SNESGetFunction()
2205 @*/
2206 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
2207 {
2208   PetscFunctionBegin;
2209   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2210   PetscValidPointer(x,2);
2211   *x = snes->vec_sol_update;
2212   PetscFunctionReturn(0);
2213 }
2214 
2215 #undef __FUNCT__
2216 #define __FUNCT__ "SNESGetFunction"
2217 /*@C
2218    SNESGetFunction - Returns the vector where the function is stored.
2219 
2220    Not Collective, but Vec is parallel if SNES is parallel
2221 
2222    Input Parameter:
2223 .  snes - the SNES context
2224 
2225    Output Parameter:
2226 +  r - the function (or PETSC_NULL)
2227 .  func - the function (or PETSC_NULL)
2228 -  ctx - the function context (or PETSC_NULL)
2229 
2230    Level: advanced
2231 
2232 .keywords: SNES, nonlinear, get, function
2233 
2234 .seealso: SNESSetFunction(), SNESGetSolution()
2235 @*/
2236 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
2237 {
2238   PetscFunctionBegin;
2239   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2240   if (r)    *r    = snes->vec_func;
2241   if (func) *func = snes->ops->computefunction;
2242   if (ctx)  *ctx  = snes->funP;
2243   PetscFunctionReturn(0);
2244 }
2245 
2246 #undef __FUNCT__
2247 #define __FUNCT__ "SNESSetOptionsPrefix"
2248 /*@C
2249    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2250    SNES options in the database.
2251 
2252    Collective on SNES
2253 
2254    Input Parameter:
2255 +  snes - the SNES context
2256 -  prefix - the prefix to prepend to all option names
2257 
2258    Notes:
2259    A hyphen (-) must NOT be given at the beginning of the prefix name.
2260    The first character of all runtime options is AUTOMATICALLY the hyphen.
2261 
2262    Level: advanced
2263 
2264 .keywords: SNES, set, options, prefix, database
2265 
2266 .seealso: SNESSetFromOptions()
2267 @*/
2268 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
2269 {
2270   PetscErrorCode ierr;
2271 
2272   PetscFunctionBegin;
2273   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2274   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2275   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2276   PetscFunctionReturn(0);
2277 }
2278 
2279 #undef __FUNCT__
2280 #define __FUNCT__ "SNESAppendOptionsPrefix"
2281 /*@C
2282    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2283    SNES options in the database.
2284 
2285    Collective on SNES
2286 
2287    Input Parameters:
2288 +  snes - the SNES context
2289 -  prefix - the prefix to prepend to all option names
2290 
2291    Notes:
2292    A hyphen (-) must NOT be given at the beginning of the prefix name.
2293    The first character of all runtime options is AUTOMATICALLY the hyphen.
2294 
2295    Level: advanced
2296 
2297 .keywords: SNES, append, options, prefix, database
2298 
2299 .seealso: SNESGetOptionsPrefix()
2300 @*/
2301 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
2302 {
2303   PetscErrorCode ierr;
2304 
2305   PetscFunctionBegin;
2306   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2307   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2308   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
2309   PetscFunctionReturn(0);
2310 }
2311 
2312 #undef __FUNCT__
2313 #define __FUNCT__ "SNESGetOptionsPrefix"
2314 /*@C
2315    SNESGetOptionsPrefix - Sets the prefix used for searching for all
2316    SNES options in the database.
2317 
2318    Not Collective
2319 
2320    Input Parameter:
2321 .  snes - the SNES context
2322 
2323    Output Parameter:
2324 .  prefix - pointer to the prefix string used
2325 
2326    Notes: On the fortran side, the user should pass in a string 'prefix' of
2327    sufficient length to hold the prefix.
2328 
2329    Level: advanced
2330 
2331 .keywords: SNES, get, options, prefix, database
2332 
2333 .seealso: SNESAppendOptionsPrefix()
2334 @*/
2335 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
2336 {
2337   PetscErrorCode ierr;
2338 
2339   PetscFunctionBegin;
2340   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2341   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
2342   PetscFunctionReturn(0);
2343 }
2344 
2345 
2346 #undef __FUNCT__
2347 #define __FUNCT__ "SNESRegister"
2348 /*@C
2349   SNESRegister - See SNESRegisterDynamic()
2350 
2351   Level: advanced
2352 @*/
2353 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2354 {
2355   char           fullname[PETSC_MAX_PATH_LEN];
2356   PetscErrorCode ierr;
2357 
2358   PetscFunctionBegin;
2359   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2360   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2361   PetscFunctionReturn(0);
2362 }
2363 
2364 #undef __FUNCT__
2365 #define __FUNCT__ "SNESTestLocalMin"
2366 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2367 {
2368   PetscErrorCode ierr;
2369   PetscInt       N,i,j;
2370   Vec            u,uh,fh;
2371   PetscScalar    value;
2372   PetscReal      norm;
2373 
2374   PetscFunctionBegin;
2375   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2376   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2377   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2378 
2379   /* currently only works for sequential */
2380   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2381   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2382   for (i=0; i<N; i++) {
2383     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2384     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2385     for (j=-10; j<11; j++) {
2386       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2387       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2388       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2389       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2390       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2391       value = -value;
2392       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2393     }
2394   }
2395   ierr = VecDestroy(uh);CHKERRQ(ierr);
2396   ierr = VecDestroy(fh);CHKERRQ(ierr);
2397   PetscFunctionReturn(0);
2398 }
2399 
2400 #undef __FUNCT__
2401 #define __FUNCT__ "SNESKSPSetUseEW"
2402 /*@
2403    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
2404    computing relative tolerance for linear solvers within an inexact
2405    Newton method.
2406 
2407    Collective on SNES
2408 
2409    Input Parameters:
2410 +  snes - SNES context
2411 -  flag - PETSC_TRUE or PETSC_FALSE
2412 
2413    Notes:
2414    Currently, the default is to use a constant relative tolerance for
2415    the inner linear solvers.  Alternatively, one can use the
2416    Eisenstat-Walker method, where the relative convergence tolerance
2417    is reset at each Newton iteration according progress of the nonlinear
2418    solver.
2419 
2420    Level: advanced
2421 
2422    Reference:
2423    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2424    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2425 
2426 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2427 
2428 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2429 @*/
2430 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag)
2431 {
2432   PetscFunctionBegin;
2433   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2434   snes->ksp_ewconv = flag;
2435   PetscFunctionReturn(0);
2436 }
2437 
2438 #undef __FUNCT__
2439 #define __FUNCT__ "SNESKSPGetUseEW"
2440 /*@
2441    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
2442    for computing relative tolerance for linear solvers within an
2443    inexact Newton method.
2444 
2445    Not Collective
2446 
2447    Input Parameter:
2448 .  snes - SNES context
2449 
2450    Output Parameter:
2451 .  flag - PETSC_TRUE or PETSC_FALSE
2452 
2453    Notes:
2454    Currently, the default is to use a constant relative tolerance for
2455    the inner linear solvers.  Alternatively, one can use the
2456    Eisenstat-Walker method, where the relative convergence tolerance
2457    is reset at each Newton iteration according progress of the nonlinear
2458    solver.
2459 
2460    Level: advanced
2461 
2462    Reference:
2463    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2464    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
2465 
2466 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
2467 
2468 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
2469 @*/
2470 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag)
2471 {
2472   PetscFunctionBegin;
2473   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2474   PetscValidPointer(flag,2);
2475   *flag = snes->ksp_ewconv;
2476   PetscFunctionReturn(0);
2477 }
2478 
2479 #undef __FUNCT__
2480 #define __FUNCT__ "SNESKSPSetParametersEW"
2481 /*@
2482    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
2483    convergence criteria for the linear solvers within an inexact
2484    Newton method.
2485 
2486    Collective on SNES
2487 
2488    Input Parameters:
2489 +    snes - SNES context
2490 .    version - version 1, 2 (default is 2) or 3
2491 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2492 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2493 .    gamma - multiplicative factor for version 2 rtol computation
2494              (0 <= gamma2 <= 1)
2495 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2496 .    alpha2 - power for safeguard
2497 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2498 
2499    Note:
2500    Version 3 was contributed by Luis Chacon, June 2006.
2501 
2502    Use PETSC_DEFAULT to retain the default for any of the parameters.
2503 
2504    Level: advanced
2505 
2506    Reference:
2507    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
2508    inexact Newton method", Utah State University Math. Stat. Dept. Res.
2509    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
2510 
2511 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
2512 
2513 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
2514 @*/
2515 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
2516 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
2517 {
2518   SNESKSPEW *kctx;
2519   PetscFunctionBegin;
2520   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2521   kctx = (SNESKSPEW*)snes->kspconvctx;
2522   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2523 
2524   if (version != PETSC_DEFAULT)   kctx->version   = version;
2525   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
2526   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
2527   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
2528   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
2529   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
2530   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
2531 
2532   if (kctx->version < 1 || kctx->version > 3) {
2533     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
2534   }
2535   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
2536     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
2537   }
2538   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
2539     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
2540   }
2541   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
2542     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
2543   }
2544   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
2545     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
2546   }
2547   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
2548     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
2549   }
2550   PetscFunctionReturn(0);
2551 }
2552 
2553 #undef __FUNCT__
2554 #define __FUNCT__ "SNESKSPGetParametersEW"
2555 /*@
2556    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
2557    convergence criteria for the linear solvers within an inexact
2558    Newton method.
2559 
2560    Not Collective
2561 
2562    Input Parameters:
2563      snes - SNES context
2564 
2565    Output Parameters:
2566 +    version - version 1, 2 (default is 2) or 3
2567 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
2568 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
2569 .    gamma - multiplicative factor for version 2 rtol computation
2570              (0 <= gamma2 <= 1)
2571 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
2572 .    alpha2 - power for safeguard
2573 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
2574 
2575    Level: advanced
2576 
2577 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
2578 
2579 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
2580 @*/
2581 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
2582 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
2583 {
2584   SNESKSPEW *kctx;
2585   PetscFunctionBegin;
2586   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2587   kctx = (SNESKSPEW*)snes->kspconvctx;
2588   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2589   if(version)   *version   = kctx->version;
2590   if(rtol_0)    *rtol_0    = kctx->rtol_0;
2591   if(rtol_max)  *rtol_max  = kctx->rtol_max;
2592   if(gamma)     *gamma     = kctx->gamma;
2593   if(alpha)     *alpha     = kctx->alpha;
2594   if(alpha2)    *alpha2    = kctx->alpha2;
2595   if(threshold) *threshold = kctx->threshold;
2596   PetscFunctionReturn(0);
2597 }
2598 
2599 #undef __FUNCT__
2600 #define __FUNCT__ "SNESKSPEW_PreSolve"
2601 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
2602 {
2603   PetscErrorCode ierr;
2604   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2605   PetscReal      rtol=PETSC_DEFAULT,stol;
2606 
2607   PetscFunctionBegin;
2608   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2609   if (!snes->iter) { /* first time in, so use the original user rtol */
2610     rtol = kctx->rtol_0;
2611   } else {
2612     if (kctx->version == 1) {
2613       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
2614       if (rtol < 0.0) rtol = -rtol;
2615       stol = pow(kctx->rtol_last,kctx->alpha2);
2616       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2617     } else if (kctx->version == 2) {
2618       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2619       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
2620       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
2621     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
2622       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
2623       /* safeguard: avoid sharp decrease of rtol */
2624       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
2625       stol = PetscMax(rtol,stol);
2626       rtol = PetscMin(kctx->rtol_0,stol);
2627       /* safeguard: avoid oversolving */
2628       stol = kctx->gamma*(snes->ttol)/snes->norm;
2629       stol = PetscMax(rtol,stol);
2630       rtol = PetscMin(kctx->rtol_0,stol);
2631     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
2632   }
2633   /* safeguard: avoid rtol greater than one */
2634   rtol = PetscMin(rtol,kctx->rtol_max);
2635   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2636   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
2637   PetscFunctionReturn(0);
2638 }
2639 
2640 #undef __FUNCT__
2641 #define __FUNCT__ "SNESKSPEW_PostSolve"
2642 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
2643 {
2644   PetscErrorCode ierr;
2645   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
2646   PCSide         pcside;
2647   Vec            lres;
2648 
2649   PetscFunctionBegin;
2650   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
2651   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
2652   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
2653   if (kctx->version == 1) {
2654     ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr);
2655     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
2656       /* KSP residual is true linear residual */
2657       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
2658     } else {
2659       /* KSP residual is preconditioned residual */
2660       /* compute true linear residual norm */
2661       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
2662       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
2663       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
2664       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
2665       ierr = VecDestroy(lres);CHKERRQ(ierr);
2666     }
2667   }
2668   PetscFunctionReturn(0);
2669 }
2670 
2671 #undef __FUNCT__
2672 #define __FUNCT__ "SNES_KSPSolve"
2673 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
2674 {
2675   PetscErrorCode ierr;
2676 
2677   PetscFunctionBegin;
2678   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
2679   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
2680   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
2681   PetscFunctionReturn(0);
2682 }
2683