xref: /petsc/src/snes/interface/snes.c (revision fde0ff241146f0ea9c891c46693ff6f0a90bbe08)
1 
2 #include <private/snesimpl.h>      /*I "petscsnes.h"  I*/
3 
4 PetscBool  SNESRegisterAllCalled = PETSC_FALSE;
5 PetscFList SNESList              = PETSC_NULL;
6 
7 /* Logging support */
8 PetscClassId  SNES_CLASSID;
9 PetscLogEvent  SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval, SNES_GSEval;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "SNESDMComputeJacobian"
13 /*
14     Translates from a SNES call to a DM call in computing a Jacobian
15 */
16 PetscErrorCode SNESDMComputeJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
17 {
18   PetscErrorCode ierr;
19   DM             dm;
20 
21   PetscFunctionBegin;
22   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
23   ierr = DMComputeJacobian(dm,X,*J,*B,flag);CHKERRQ(ierr);
24   PetscFunctionReturn(0);
25 }
26 
27 #undef __FUNCT__
28 #define __FUNCT__ "SNESSetErrorIfNotConverged"
29 /*@
30    SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged.
31 
32    Logically Collective on SNES
33 
34    Input Parameters:
35 +  snes - iterative context obtained from SNESCreate()
36 -  flg - PETSC_TRUE indicates you want the error generated
37 
38    Options database keys:
39 .  -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
40 
41    Level: intermediate
42 
43    Notes:
44     Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve()
45     to determine if it has converged.
46 
47 .keywords: SNES, set, initial guess, nonzero
48 
49 .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
50 @*/
51 PetscErrorCode  SNESSetErrorIfNotConverged(SNES snes,PetscBool  flg)
52 {
53   PetscFunctionBegin;
54   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
55   PetscValidLogicalCollectiveBool(snes,flg,2);
56   snes->errorifnotconverged = flg;
57 
58   PetscFunctionReturn(0);
59 }
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "SNESGetErrorIfNotConverged"
63 /*@
64    SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge?
65 
66    Not Collective
67 
68    Input Parameter:
69 .  snes - iterative context obtained from SNESCreate()
70 
71    Output Parameter:
72 .  flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
73 
74    Level: intermediate
75 
76 .keywords: SNES, set, initial guess, nonzero
77 
78 .seealso:  SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
79 @*/
80 PetscErrorCode  SNESGetErrorIfNotConverged(SNES snes,PetscBool  *flag)
81 {
82   PetscFunctionBegin;
83   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
84   PetscValidPointer(flag,2);
85   *flag = snes->errorifnotconverged;
86   PetscFunctionReturn(0);
87 }
88 
89 #undef __FUNCT__
90 #define __FUNCT__ "SNESSetFunctionDomainError"
91 /*@
92    SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not
93      in the functions domain. For example, negative pressure.
94 
95    Logically Collective on SNES
96 
97    Input Parameters:
98 .  SNES - the SNES context
99 
100    Level: advanced
101 
102 .keywords: SNES, view
103 
104 .seealso: SNESCreate(), SNESSetFunction()
105 @*/
106 PetscErrorCode  SNESSetFunctionDomainError(SNES snes)
107 {
108   PetscFunctionBegin;
109   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
110   snes->domainerror = PETSC_TRUE;
111   PetscFunctionReturn(0);
112 }
113 
114 #undef __FUNCT__
115 #define __FUNCT__ "SNESView"
116 /*@C
117    SNESView - Prints the SNES data structure.
118 
119    Collective on SNES
120 
121    Input Parameters:
122 +  SNES - the SNES context
123 -  viewer - visualization context
124 
125    Options Database Key:
126 .  -snes_view - Calls SNESView() at end of SNESSolve()
127 
128    Notes:
129    The available visualization contexts include
130 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
131 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
132          output where only the first processor opens
133          the file.  All other processors send their
134          data to the first processor to print.
135 
136    The user can open an alternative visualization context with
137    PetscViewerASCIIOpen() - output to a specified file.
138 
139    Level: beginner
140 
141 .keywords: SNES, view
142 
143 .seealso: PetscViewerASCIIOpen()
144 @*/
145 PetscErrorCode  SNESView(SNES snes,PetscViewer viewer)
146 {
147   SNESKSPEW           *kctx;
148   PetscErrorCode      ierr;
149   KSP                 ksp;
150   PetscBool           iascii,isstring;
151 
152   PetscFunctionBegin;
153   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
154   if (!viewer) {
155     ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
156   }
157   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
158   PetscCheckSameComm(snes,1,viewer,2);
159 
160   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
161   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
162   if (iascii) {
163     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer,"SNES Object");CHKERRQ(ierr);
164     if (snes->ops->view) {
165       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
166       ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);
167       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
168     }
169     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
170     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
171                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
172     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
173     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
174     if (snes->ksp_ewconv) {
175       kctx = (SNESKSPEW *)snes->kspconvctx;
176       if (kctx) {
177         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
178         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
179         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
180       }
181     }
182     if (snes->lagpreconditioner == -1) {
183       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioned is never rebuilt\n");CHKERRQ(ierr);
184     } else if (snes->lagpreconditioner > 1) {
185       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);CHKERRQ(ierr);
186     }
187     if (snes->lagjacobian == -1) {
188       ierr = PetscViewerASCIIPrintf(viewer,"  Jacobian is never rebuilt\n");CHKERRQ(ierr);
189     } else if (snes->lagjacobian > 1) {
190       ierr = PetscViewerASCIIPrintf(viewer,"  Jacobian is rebuilt every %D SNES iterations\n",snes->lagjacobian);CHKERRQ(ierr);
191     }
192   } else if (isstring) {
193     const char *type;
194     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
195     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
196   }
197   if (snes->pc && snes->usespc) {
198     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
199     ierr = SNESView(snes->pc, viewer);CHKERRQ(ierr);
200     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
201   }
202   if (snes->usesksp) {
203     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
204     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
205     ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
206     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
207   }
208   PetscFunctionReturn(0);
209 }
210 
211 /*
212   We retain a list of functions that also take SNES command
213   line options. These are called at the end SNESSetFromOptions()
214 */
215 #define MAXSETFROMOPTIONS 5
216 static PetscInt numberofsetfromoptions;
217 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "SNESAddOptionsChecker"
221 /*@C
222   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
223 
224   Not Collective
225 
226   Input Parameter:
227 . snescheck - function that checks for options
228 
229   Level: developer
230 
231 .seealso: SNESSetFromOptions()
232 @*/
233 PetscErrorCode  SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
234 {
235   PetscFunctionBegin;
236   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
237     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
238   }
239   othersetfromoptions[numberofsetfromoptions++] = snescheck;
240   PetscFunctionReturn(0);
241 }
242 
243 extern PetscErrorCode  SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
244 
245 #undef __FUNCT__
246 #define __FUNCT__ "SNESSetUpMatrixFree_Private"
247 static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool  hasOperator, PetscInt version)
248 {
249   Mat            J;
250   KSP            ksp;
251   PC             pc;
252   PetscBool      match;
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
257 
258   if(!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
259     Mat A = snes->jacobian, B = snes->jacobian_pre;
260     ierr = MatGetVecs(A ? A : B, PETSC_NULL,&snes->vec_func);CHKERRQ(ierr);
261   }
262 
263   if (version == 1) {
264     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
265     ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr);
266     ierr = MatSetFromOptions(J);CHKERRQ(ierr);
267   } else if (version == 2) {
268     if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
269 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
270     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);CHKERRQ(ierr);
271 #else
272     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
273 #endif
274   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");
275 
276   ierr = PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);CHKERRQ(ierr);
277   if (hasOperator) {
278     /* This version replaces the user provided Jacobian matrix with a
279        matrix-free version but still employs the user-provided preconditioner matrix. */
280     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
281   } else {
282     /* This version replaces both the user-provided Jacobian and the user-
283        provided preconditioner matrix with the default matrix free version. */
284     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
285     /* Force no preconditioner */
286     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
287     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
288     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&match);CHKERRQ(ierr);
289     if (!match) {
290       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
291       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
292     }
293   }
294   ierr = MatDestroy(&J);CHKERRQ(ierr);
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "SNESSetFromOptions"
300 /*@
301    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
302 
303    Collective on SNES
304 
305    Input Parameter:
306 .  snes - the SNES context
307 
308    Options Database Keys:
309 +  -snes_type <type> - ls, tr, ngmres, ncg, richardson, qn, vi, fas
310 .  -snes_stol - convergence tolerance in terms of the norm
311                 of the change in the solution between steps
312 .  -snes_atol <abstol> - absolute tolerance of residual norm
313 .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
314 .  -snes_max_it <max_it> - maximum number of iterations
315 .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
316 .  -snes_max_fail <max_fail> - maximum number of line search failures allowed before stopping, default is none
317 .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
318 .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
319 .  -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
320 .  -snes_trtol <trtol> - trust region tolerance
321 .  -snes_no_convergence_test - skip convergence test in nonlinear
322                                solver; hence iterations will continue until max_it
323                                or some other criterion is reached. Saves expense
324                                of convergence test
325 .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
326                                        filename given prints to stdout
327 .  -snes_monitor_solution - plots solution at each iteration
328 .  -snes_monitor_residual - plots residual (not its norm) at each iteration
329 .  -snes_monitor_solution_update - plots update to solution at each iteration
330 .  -snes_monitor_draw - plots residual norm at each iteration
331 .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
332 .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
333 -  -snes_converged_reason - print the reason for convergence/divergence after each solve
334 
335     Options Database for Eisenstat-Walker method:
336 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
337 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
338 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
339 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
340 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
341 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
342 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
343 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
344 
345    Notes:
346    To see all options, run your program with the -help option or consult
347    the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A>.
348 
349    Level: beginner
350 
351 .keywords: SNES, nonlinear, set, options, database
352 
353 .seealso: SNESSetOptionsPrefix()
354 @*/
355 PetscErrorCode  SNESSetFromOptions(SNES snes)
356 {
357   PetscBool               flg,set,mf,mf_operator,pcset;
358   PetscInt                i,indx,lag,mf_version,grids;
359   MatStructure            matflag;
360   const char              *deft = SNESLS;
361   const char              *convtests[] = {"default","skip"};
362   SNESKSPEW               *kctx = NULL;
363   char                    type[256], monfilename[PETSC_MAX_PATH_LEN];
364   const char              *optionsprefix;
365   PetscViewer             monviewer;
366   PetscErrorCode          ierr;
367 
368   PetscFunctionBegin;
369   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
370 
371   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
372   ierr = PetscObjectOptionsBegin((PetscObject)snes);CHKERRQ(ierr);
373     if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; }
374     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
375     if (flg) {
376       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
377     } else if (!((PetscObject)snes)->type_name) {
378       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
379     }
380     /* not used here, but called so will go into help messaage */
381     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
382 
383     ierr = PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
384     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
385 
386     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
387     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
388     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
389     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
390     ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr);
391     ierr = PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,PETSC_NULL);CHKERRQ(ierr);
392 
393     ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr);
394     if (flg) {
395       ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr);
396     }
397     ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr);
398     if (flg) {
399       ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr);
400     }
401     ierr = PetscOptionsInt("-snes_grid_sequence","Use grid sequencing to generate initial guess","SNESSetGridSequence",snes->gridsequence,&grids,&flg);CHKERRQ(ierr);
402     if (flg) {
403       ierr = SNESSetGridSequence(snes,grids);CHKERRQ(ierr);
404     }
405 
406     ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr);
407     if (flg) {
408       switch (indx) {
409       case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break;
410       case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);    break;
411       }
412     }
413 
414     ierr = PetscOptionsBool("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",snes->printreason,&snes->printreason,PETSC_NULL);CHKERRQ(ierr);
415 
416     kctx = (SNESKSPEW *)snes->kspconvctx;
417 
418     ierr = PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr);
419 
420     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
421     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
422     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
423     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
424     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
425     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
426     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
427 
428     flg  = PETSC_FALSE;
429     ierr = PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
430     if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);}
431 
432     ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
433     if (flg) {
434       ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,monfilename,&monviewer);CHKERRQ(ierr);
435       ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
436     }
437 
438     ierr = PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
439     if (flg) {
440       ierr = SNESMonitorSet(snes,SNESMonitorRange,0,0);CHKERRQ(ierr);
441     }
442 
443     ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
444     if (flg) {
445       ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,monfilename,&monviewer);CHKERRQ(ierr);
446       ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr);
447     }
448 
449     ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
450     if (flg) {
451       ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,monfilename,&monviewer);CHKERRQ(ierr);
452       ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
453     }
454 
455     ierr = PetscOptionsString("-snes_monitor_python","Use Python function","SNESMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
456     if (flg) {ierr = PetscPythonMonitorSet((PetscObject)snes,monfilename);CHKERRQ(ierr);}
457 
458     flg  = PETSC_FALSE;
459     ierr = PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
460     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);}
461     flg  = PETSC_FALSE;
462     ierr = PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
463     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);}
464     flg  = PETSC_FALSE;
465     ierr = PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
466     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);}
467     flg  = PETSC_FALSE;
468     ierr = PetscOptionsBool("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
469     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
470     flg  = PETSC_FALSE;
471     ierr = PetscOptionsBool("-snes_monitor_range_draw","Plot function range at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
472     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLGRange,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
473 
474     flg  = PETSC_FALSE;
475     ierr = PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
476     if (flg) {
477       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
478       ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr);
479     }
480 
481     mf = mf_operator = PETSC_FALSE;
482     flg = PETSC_FALSE;
483     ierr = PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf_operator,&flg);CHKERRQ(ierr);
484     if (flg && mf_operator) {
485       snes->mf_operator = PETSC_TRUE;
486       mf = PETSC_TRUE;
487     }
488     flg = PETSC_FALSE;
489     ierr = PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf,&flg);CHKERRQ(ierr);
490     if (!flg && mf_operator) mf = PETSC_TRUE;
491     mf_version = 1;
492     ierr = PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",mf_version,&mf_version,0);CHKERRQ(ierr);
493 
494 
495     /* GS Options */
496     ierr = PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",snes->gssweeps,&snes->gssweeps,PETSC_NULL);CHKERRQ(ierr);
497 
498     /* line search options */
499     ierr = PetscOptionsReal("-snes_ls_alpha","Constant function norm must decrease by","None",snes->ls_alpha,&snes->ls_alpha,0);CHKERRQ(ierr);
500     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",snes->maxstep,&snes->maxstep,0);CHKERRQ(ierr);
501     ierr = PetscOptionsReal("-snes_ls_steptol","Minimum lambda allowed","None",snes->steptol,&snes->steptol,0);CHKERRQ(ierr);
502     ierr = PetscOptionsReal("-snes_ls_damping","Damping parameter","SNES",snes->damping,&snes->damping,&flg);CHKERRQ(ierr);
503     ierr = PetscOptionsInt("-snes_ls_it"      ,"Line search iterations","SNES",snes->ls_its,&snes->ls_its,&flg);CHKERRQ(ierr);
504     ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",snes->ls_monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
505     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
506     ierr = PetscOptionsEnum("-snes_ls","Line search used","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)snes->ls_type,(PetscEnum*)&indx,&flg);CHKERRQ(ierr);
507     if (flg) {
508       ierr = SNESLineSearchSetType(snes,(SNESLineSearchType)indx);CHKERRQ(ierr);
509     }
510     flg = snes->ops->precheckstep == SNESLineSearchPreCheckPicard ? PETSC_TRUE : PETSC_FALSE;
511     ierr = PetscOptionsBool("-snes_ls_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr);
512     if (set) {
513       if (flg) {
514         snes->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */
515         ierr = PetscOptionsReal("-snes_ls_precheck_picard_angle","Maximum angle at which to activate the correction","none",snes->precheck_picard_angle,&snes->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr);
516         ierr = SNESLineSearchSetPreCheck(snes,SNESLineSearchPreCheckPicard,&snes->precheck_picard_angle);CHKERRQ(ierr);
517       } else {
518         ierr = SNESLineSearchSetPreCheck(snes,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
519       }
520     }
521 
522     for(i = 0; i < numberofsetfromoptions; i++) {
523       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
524     }
525 
526     if (snes->ops->setfromoptions) {
527       ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr);
528     }
529 
530     /* process any options handlers added with PetscObjectAddOptionsHandler() */
531     ierr = PetscObjectProcessOptionsHandlers((PetscObject)snes);CHKERRQ(ierr);
532   ierr = PetscOptionsEnd();CHKERRQ(ierr);
533 
534   if (mf) { ierr = SNESSetUpMatrixFree_Private(snes, mf_operator, mf_version);CHKERRQ(ierr); }
535 
536   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
537   ierr = KSPGetOperators(snes->ksp,PETSC_NULL,PETSC_NULL,&matflag);
538   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,matflag);CHKERRQ(ierr);
539   ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr);
540 
541   /* if someone has set the SNES PC type, create it. */
542   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
543   ierr = PetscOptionsHasName(optionsprefix, "-npc_snes_type", &pcset);CHKERRQ(ierr);
544   if (pcset && (!snes->pc)) {
545     ierr = SNESGetPC(snes, &snes->pc);CHKERRQ(ierr);
546   }
547   if (snes->pc) {
548     ierr = SNESSetOptionsPrefix(snes->pc, optionsprefix);CHKERRQ(ierr);
549     ierr = SNESAppendOptionsPrefix(snes->pc, "npc_");CHKERRQ(ierr);
550     ierr = SNESSetDM(snes->pc, snes->dm);CHKERRQ(ierr);
551     ierr = SNESSetGS(snes->pc, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
552     /* Should we make a duplicate vector and matrix? Leave the DM to make it? */
553     ierr = SNESSetFunction(snes->pc, snes->vec_func, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
554     ierr = SNESSetJacobian(snes->pc, snes->jacobian, snes->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
555     ierr = SNESSetFromOptions(snes->pc);CHKERRQ(ierr);
556   }
557   PetscFunctionReturn(0);
558 }
559 
560 #undef __FUNCT__
561 #define __FUNCT__ "SNESSetComputeApplicationContext"
562 /*@
563    SNESSetComputeApplicationContext - Sets an optional function to compute a user-defined context for
564    the nonlinear solvers.
565 
566    Logically Collective on SNES
567 
568    Input Parameters:
569 +  snes - the SNES context
570 .  compute - function to compute the context
571 -  destroy - function to destroy the context
572 
573    Level: intermediate
574 
575 .keywords: SNES, nonlinear, set, application, context
576 
577 .seealso: SNESGetApplicationContext(), SNESSetComputeApplicationContext(), SNESGetApplicationContext()
578 @*/
579 PetscErrorCode  SNESSetComputeApplicationContext(SNES snes,PetscErrorCode (*compute)(SNES,void**),PetscErrorCode (*destroy)(void**))
580 {
581   PetscFunctionBegin;
582   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
583   snes->ops->usercompute = compute;
584   snes->ops->userdestroy = destroy;
585   PetscFunctionReturn(0);
586 }
587 
588 #undef __FUNCT__
589 #define __FUNCT__ "SNESSetApplicationContext"
590 /*@
591    SNESSetApplicationContext - Sets the optional user-defined context for
592    the nonlinear solvers.
593 
594    Logically Collective on SNES
595 
596    Input Parameters:
597 +  snes - the SNES context
598 -  usrP - optional user context
599 
600    Level: intermediate
601 
602 .keywords: SNES, nonlinear, set, application, context
603 
604 .seealso: SNESGetApplicationContext(), SNESSetApplicationContext()
605 @*/
606 PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
607 {
608   PetscErrorCode ierr;
609   KSP            ksp;
610 
611   PetscFunctionBegin;
612   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
613   ierr       = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
614   ierr       = KSPSetApplicationContext(ksp,usrP);CHKERRQ(ierr);
615   snes->user = usrP;
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "SNESGetApplicationContext"
621 /*@
622    SNESGetApplicationContext - Gets the user-defined context for the
623    nonlinear solvers.
624 
625    Not Collective
626 
627    Input Parameter:
628 .  snes - SNES context
629 
630    Output Parameter:
631 .  usrP - user context
632 
633    Level: intermediate
634 
635 .keywords: SNES, nonlinear, get, application, context
636 
637 .seealso: SNESSetApplicationContext()
638 @*/
639 PetscErrorCode  SNESGetApplicationContext(SNES snes,void *usrP)
640 {
641   PetscFunctionBegin;
642   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
643   *(void**)usrP = snes->user;
644   PetscFunctionReturn(0);
645 }
646 
647 #undef __FUNCT__
648 #define __FUNCT__ "SNESGetIterationNumber"
649 /*@
650    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
651    at this time.
652 
653    Not Collective
654 
655    Input Parameter:
656 .  snes - SNES context
657 
658    Output Parameter:
659 .  iter - iteration number
660 
661    Notes:
662    For example, during the computation of iteration 2 this would return 1.
663 
664    This is useful for using lagged Jacobians (where one does not recompute the
665    Jacobian at each SNES iteration). For example, the code
666 .vb
667       ierr = SNESGetIterationNumber(snes,&it);
668       if (!(it % 2)) {
669         [compute Jacobian here]
670       }
671 .ve
672    can be used in your ComputeJacobian() function to cause the Jacobian to be
673    recomputed every second SNES iteration.
674 
675    Level: intermediate
676 
677 .keywords: SNES, nonlinear, get, iteration, number,
678 
679 .seealso:   SNESGetFunctionNorm(), SNESGetLinearSolveIterations()
680 @*/
681 PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt* iter)
682 {
683   PetscFunctionBegin;
684   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
685   PetscValidIntPointer(iter,2);
686   *iter = snes->iter;
687   PetscFunctionReturn(0);
688 }
689 
690 #undef __FUNCT__
691 #define __FUNCT__ "SNESGetFunctionNorm"
692 /*@
693    SNESGetFunctionNorm - Gets the norm of the current function that was set
694    with SNESSSetFunction().
695 
696    Collective on SNES
697 
698    Input Parameter:
699 .  snes - SNES context
700 
701    Output Parameter:
702 .  fnorm - 2-norm of function
703 
704    Level: intermediate
705 
706 .keywords: SNES, nonlinear, get, function, norm
707 
708 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations()
709 @*/
710 PetscErrorCode  SNESGetFunctionNorm(SNES snes,PetscReal *fnorm)
711 {
712   PetscFunctionBegin;
713   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
714   PetscValidScalarPointer(fnorm,2);
715   *fnorm = snes->norm;
716   PetscFunctionReturn(0);
717 }
718 
719 #undef __FUNCT__
720 #define __FUNCT__ "SNESGetNonlinearStepFailures"
721 /*@
722    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
723    attempted by the nonlinear solver.
724 
725    Not Collective
726 
727    Input Parameter:
728 .  snes - SNES context
729 
730    Output Parameter:
731 .  nfails - number of unsuccessful steps attempted
732 
733    Notes:
734    This counter is reset to zero for each successive call to SNESSolve().
735 
736    Level: intermediate
737 
738 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
739 
740 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
741           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
742 @*/
743 PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails)
744 {
745   PetscFunctionBegin;
746   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
747   PetscValidIntPointer(nfails,2);
748   *nfails = snes->numFailures;
749   PetscFunctionReturn(0);
750 }
751 
752 #undef __FUNCT__
753 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures"
754 /*@
755    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
756    attempted by the nonlinear solver before it gives up.
757 
758    Not Collective
759 
760    Input Parameters:
761 +  snes     - SNES context
762 -  maxFails - maximum of unsuccessful steps
763 
764    Level: intermediate
765 
766 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
767 
768 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
769           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
770 @*/
771 PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
772 {
773   PetscFunctionBegin;
774   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
775   snes->maxFailures = maxFails;
776   PetscFunctionReturn(0);
777 }
778 
779 #undef __FUNCT__
780 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures"
781 /*@
782    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
783    attempted by the nonlinear solver before it gives up.
784 
785    Not Collective
786 
787    Input Parameter:
788 .  snes     - SNES context
789 
790    Output Parameter:
791 .  maxFails - maximum of unsuccessful steps
792 
793    Level: intermediate
794 
795 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
796 
797 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
798           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
799 
800 @*/
801 PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
802 {
803   PetscFunctionBegin;
804   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
805   PetscValidIntPointer(maxFails,2);
806   *maxFails = snes->maxFailures;
807   PetscFunctionReturn(0);
808 }
809 
810 #undef __FUNCT__
811 #define __FUNCT__ "SNESGetNumberFunctionEvals"
812 /*@
813    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
814      done by SNES.
815 
816    Not Collective
817 
818    Input Parameter:
819 .  snes     - SNES context
820 
821    Output Parameter:
822 .  nfuncs - number of evaluations
823 
824    Level: intermediate
825 
826 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
827 
828 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures()
829 @*/
830 PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
831 {
832   PetscFunctionBegin;
833   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
834   PetscValidIntPointer(nfuncs,2);
835   *nfuncs = snes->nfuncs;
836   PetscFunctionReturn(0);
837 }
838 
839 #undef __FUNCT__
840 #define __FUNCT__ "SNESGetLinearSolveFailures"
841 /*@
842    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
843    linear solvers.
844 
845    Not Collective
846 
847    Input Parameter:
848 .  snes - SNES context
849 
850    Output Parameter:
851 .  nfails - number of failed solves
852 
853    Notes:
854    This counter is reset to zero for each successive call to SNESSolve().
855 
856    Level: intermediate
857 
858 .keywords: SNES, nonlinear, get, number, unsuccessful, steps
859 
860 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
861 @*/
862 PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails)
863 {
864   PetscFunctionBegin;
865   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
866   PetscValidIntPointer(nfails,2);
867   *nfails = snes->numLinearSolveFailures;
868   PetscFunctionReturn(0);
869 }
870 
871 #undef __FUNCT__
872 #define __FUNCT__ "SNESSetMaxLinearSolveFailures"
873 /*@
874    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
875    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
876 
877    Logically Collective on SNES
878 
879    Input Parameters:
880 +  snes     - SNES context
881 -  maxFails - maximum allowed linear solve failures
882 
883    Level: intermediate
884 
885    Notes: By default this is 0; that is SNES returns on the first failed linear solve
886 
887 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
888 
889 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
890 @*/
891 PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
892 {
893   PetscFunctionBegin;
894   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
895   PetscValidLogicalCollectiveInt(snes,maxFails,2);
896   snes->maxLinearSolveFailures = maxFails;
897   PetscFunctionReturn(0);
898 }
899 
900 #undef __FUNCT__
901 #define __FUNCT__ "SNESGetMaxLinearSolveFailures"
902 /*@
903    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
904      are allowed before SNES terminates
905 
906    Not Collective
907 
908    Input Parameter:
909 .  snes     - SNES context
910 
911    Output Parameter:
912 .  maxFails - maximum of unsuccessful solves allowed
913 
914    Level: intermediate
915 
916    Notes: By default this is 1; that is SNES returns on the first failed linear solve
917 
918 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
919 
920 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
921 @*/
922 PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
923 {
924   PetscFunctionBegin;
925   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
926   PetscValidIntPointer(maxFails,2);
927   *maxFails = snes->maxLinearSolveFailures;
928   PetscFunctionReturn(0);
929 }
930 
931 #undef __FUNCT__
932 #define __FUNCT__ "SNESGetLinearSolveIterations"
933 /*@
934    SNESGetLinearSolveIterations - Gets the total number of linear iterations
935    used by the nonlinear solver.
936 
937    Not Collective
938 
939    Input Parameter:
940 .  snes - SNES context
941 
942    Output Parameter:
943 .  lits - number of linear iterations
944 
945    Notes:
946    This counter is reset to zero for each successive call to SNESSolve().
947 
948    Level: intermediate
949 
950 .keywords: SNES, nonlinear, get, number, linear, iterations
951 
952 .seealso:  SNESGetIterationNumber(), SNESGetFunctionNorm(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures()
953 @*/
954 PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt* lits)
955 {
956   PetscFunctionBegin;
957   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
958   PetscValidIntPointer(lits,2);
959   *lits = snes->linear_its;
960   PetscFunctionReturn(0);
961 }
962 
963 #undef __FUNCT__
964 #define __FUNCT__ "SNESGetKSP"
965 /*@
966    SNESGetKSP - Returns the KSP context for a SNES solver.
967 
968    Not Collective, but if SNES object is parallel, then KSP object is parallel
969 
970    Input Parameter:
971 .  snes - the SNES context
972 
973    Output Parameter:
974 .  ksp - the KSP context
975 
976    Notes:
977    The user can then directly manipulate the KSP context to set various
978    options, etc.  Likewise, the user can then extract and manipulate the
979    PC contexts as well.
980 
981    Level: beginner
982 
983 .keywords: SNES, nonlinear, get, KSP, context
984 
985 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
986 @*/
987 PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
988 {
989   PetscErrorCode ierr;
990 
991   PetscFunctionBegin;
992   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
993   PetscValidPointer(ksp,2);
994 
995   if (!snes->ksp) {
996     ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr);
997     ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr);
998     ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
999   }
1000   *ksp = snes->ksp;
1001   PetscFunctionReturn(0);
1002 }
1003 
1004 #undef __FUNCT__
1005 #define __FUNCT__ "SNESSetKSP"
1006 /*@
1007    SNESSetKSP - Sets a KSP context for the SNES object to use
1008 
1009    Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
1010 
1011    Input Parameters:
1012 +  snes - the SNES context
1013 -  ksp - the KSP context
1014 
1015    Notes:
1016    The SNES object already has its KSP object, you can obtain with SNESGetKSP()
1017    so this routine is rarely needed.
1018 
1019    The KSP object that is already in the SNES object has its reference count
1020    decreased by one.
1021 
1022    Level: developer
1023 
1024 .keywords: SNES, nonlinear, get, KSP, context
1025 
1026 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
1027 @*/
1028 PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
1029 {
1030   PetscErrorCode ierr;
1031 
1032   PetscFunctionBegin;
1033   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1034   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
1035   PetscCheckSameComm(snes,1,ksp,2);
1036   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
1037   if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);}
1038   snes->ksp = ksp;
1039   PetscFunctionReturn(0);
1040 }
1041 
1042 #if 0
1043 #undef __FUNCT__
1044 #define __FUNCT__ "SNESPublish_Petsc"
1045 static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
1046 {
1047   PetscFunctionBegin;
1048   PetscFunctionReturn(0);
1049 }
1050 #endif
1051 
1052 /* -----------------------------------------------------------*/
1053 #undef __FUNCT__
1054 #define __FUNCT__ "SNESCreate"
1055 /*@
1056    SNESCreate - Creates a nonlinear solver context.
1057 
1058    Collective on MPI_Comm
1059 
1060    Input Parameters:
1061 .  comm - MPI communicator
1062 
1063    Output Parameter:
1064 .  outsnes - the new SNES context
1065 
1066    Options Database Keys:
1067 +   -snes_mf - Activates default matrix-free Jacobian-vector products,
1068                and no preconditioning matrix
1069 .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
1070                products, and a user-provided preconditioning matrix
1071                as set by SNESSetJacobian()
1072 -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
1073 
1074    Level: beginner
1075 
1076 .keywords: SNES, nonlinear, create, context
1077 
1078 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
1079 
1080 @*/
1081 PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
1082 {
1083   PetscErrorCode      ierr;
1084   SNES                snes;
1085   SNESKSPEW           *kctx;
1086 
1087   PetscFunctionBegin;
1088   PetscValidPointer(outsnes,2);
1089   *outsnes = PETSC_NULL;
1090 #ifndef PETSC_USE_DYNAMIC_LIBRARIES
1091   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
1092 #endif
1093 
1094   ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_CLASSID,0,"SNES","Nonlinear solver","SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
1095 
1096   snes->ops->converged    = SNESDefaultConverged;
1097   snes->usesksp           = PETSC_TRUE;
1098   snes->usegs             = PETSC_FALSE;
1099   snes->max_its           = 50;
1100   snes->max_funcs	  = 10000;
1101   snes->norm		  = 0.0;
1102   snes->rtol		  = 1.e-8;
1103   snes->ttol              = 0.0;
1104   snes->abstol		  = 1.e-50;
1105   snes->xtol		  = 1.e-8;
1106   snes->deltatol	  = 1.e-12;
1107   snes->nfuncs            = 0;
1108   snes->numFailures       = 0;
1109   snes->maxFailures       = 1;
1110   snes->linear_its        = 0;
1111   snes->lagjacobian       = 1;
1112   snes->lagpreconditioner = 1;
1113   snes->numbermonitors    = 0;
1114   snes->data              = 0;
1115   snes->setupcalled       = PETSC_FALSE;
1116   snes->ksp_ewconv        = PETSC_FALSE;
1117   snes->nwork             = 0;
1118   snes->work              = 0;
1119   snes->nvwork            = 0;
1120   snes->vwork             = 0;
1121   snes->conv_hist_len     = 0;
1122   snes->conv_hist_max     = 0;
1123   snes->conv_hist         = PETSC_NULL;
1124   snes->conv_hist_its     = PETSC_NULL;
1125   snes->conv_hist_reset   = PETSC_TRUE;
1126   snes->reason            = SNES_CONVERGED_ITERATING;
1127   snes->gssweeps          = 1;
1128 
1129   /* initialize the line search options */
1130   snes->ls_type           = SNES_LS_BASIC;
1131   snes->ls_its            = 1;
1132   snes->damping           = 1.0;
1133   snes->maxstep           = 1e8;
1134   snes->steptol           = 1e-12;
1135   snes->ls_alpha          = 1e-4;
1136   snes->ls_monitor        = PETSC_NULL;
1137 
1138   snes->ops->linesearch   = PETSC_NULL;
1139   snes->precheck          = PETSC_NULL;
1140   snes->ops->precheckstep = PETSC_NULL;
1141   snes->postcheck         = PETSC_NULL;
1142   snes->ops->postcheckstep= PETSC_NULL;
1143 
1144   snes->numLinearSolveFailures = 0;
1145   snes->maxLinearSolveFailures = 1;
1146 
1147   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
1148   ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr);
1149   snes->kspconvctx  = (void*)kctx;
1150   kctx->version     = 2;
1151   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
1152                              this was too large for some test cases */
1153   kctx->rtol_last   = 0.0;
1154   kctx->rtol_max    = .9;
1155   kctx->gamma       = 1.0;
1156   kctx->alpha       = .5*(1.0 + PetscSqrtReal(5.0));
1157   kctx->alpha2      = kctx->alpha;
1158   kctx->threshold   = .1;
1159   kctx->lresid_last = 0.0;
1160   kctx->norm_last   = 0.0;
1161 
1162   *outsnes = snes;
1163   PetscFunctionReturn(0);
1164 }
1165 
1166 #undef __FUNCT__
1167 #define __FUNCT__ "SNESSetFunction"
1168 /*@C
1169    SNESSetFunction - Sets the function evaluation routine and function
1170    vector for use by the SNES routines in solving systems of nonlinear
1171    equations.
1172 
1173    Logically Collective on SNES
1174 
1175    Input Parameters:
1176 +  snes - the SNES context
1177 .  r - vector to store function value
1178 .  func - function evaluation routine
1179 -  ctx - [optional] user-defined context for private data for the
1180          function evaluation routine (may be PETSC_NULL)
1181 
1182    Calling sequence of func:
1183 $    func (SNES snes,Vec x,Vec f,void *ctx);
1184 
1185 .  f - function vector
1186 -  ctx - optional user-defined function context
1187 
1188    Notes:
1189    The Newton-like methods typically solve linear systems of the form
1190 $      f'(x) x = -f(x),
1191    where f'(x) denotes the Jacobian matrix and f(x) is the function.
1192 
1193    Level: beginner
1194 
1195 .keywords: SNES, nonlinear, set, function
1196 
1197 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard()
1198 @*/
1199 PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
1200 {
1201   PetscErrorCode ierr;
1202   PetscFunctionBegin;
1203   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1204   if (r) {
1205     PetscValidHeaderSpecific(r,VEC_CLASSID,2);
1206     PetscCheckSameComm(snes,1,r,2);
1207     ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr);
1208     ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr);
1209     snes->vec_func = r;
1210   } else if (!snes->vec_func && snes->dm) {
1211     ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr);
1212   }
1213   if (func) snes->ops->computefunction = func;
1214   if (ctx)  snes->funP                 = ctx;
1215   PetscFunctionReturn(0);
1216 }
1217 
1218 
1219 #undef __FUNCT__
1220 #define __FUNCT__ "SNESSetGS"
1221 /*@C
1222    SNESSetGS - Sets the user nonlinear Gauss-Seidel routine for
1223    use with composed nonlinear solvers.
1224 
1225    Input Parameters:
1226 +  snes   - the SNES context
1227 .  gsfunc - function evaluation routine
1228 -  ctx    - [optional] user-defined context for private data for the
1229             smoother evaluation routine (may be PETSC_NULL)
1230 
1231    Calling sequence of func:
1232 $    func (SNES snes,Vec x,Vec b,void *ctx);
1233 
1234 +  X   - solution vector
1235 .  B   - RHS vector
1236 -  ctx - optional user-defined Gauss-Seidel context
1237 
1238    Notes:
1239    The GS routines are used by the composed nonlinear solver to generate
1240     a problem appropriate update to the solution, particularly FAS.
1241 
1242    Level: intermediate
1243 
1244 .keywords: SNES, nonlinear, set, Gauss-Seidel
1245 
1246 .seealso: SNESGetFunction(), SNESComputeGS()
1247 @*/
1248 PetscErrorCode SNESSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx) {
1249   PetscFunctionBegin;
1250   if (gsfunc) snes->ops->computegs = gsfunc;
1251   if (ctx) snes->gsP = ctx;
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 #undef __FUNCT__
1256 #define __FUNCT__ "SNESSetUseGS"
1257 /*@
1258    SNESSetUseGS - Toggles use of the user-set nonlinear GS solver in
1259    SNES solvers that may use a preconditioning routine.
1260 
1261    Input Parameters:
1262 +  snes   - the SNES context
1263 -  usegs  - whether to use the nonlinear GS routine or not.
1264 
1265    Notes:
1266    The behavior of the nonlinear GS solvers is that if useGS is set, then
1267    the update is constructed with the user-defined nonlinear GS method.
1268    Otherwise, the update is either formed by the user-customized
1269    preconditioner SNES, or by nonlinear richardson if both of these
1270    are not provided.
1271 
1272    Level: intermediate
1273 
1274 .keywords: SNES, nonlinear, set, Gauss-Siedel
1275 
1276 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGetUseGS()
1277 @*/
1278 PetscErrorCode SNESSetUseGS(SNES snes, PetscBool usegs) {
1279   PetscFunctionBegin;
1280   snes->usegs = usegs;
1281   PetscFunctionReturn(0);
1282 }
1283 
1284 
1285 #undef __FUNCT__
1286 #define __FUNCT__ "SNESGetUseGS"
1287 /*@
1288    SNESGetUseGS - Returns whether or not the SNES is using a
1289    provided nonlinear GS routine.
1290 
1291    Input Parameters:
1292 +  snes   - the SNES context
1293 -  usegs  - flag indicating whether or not GS is being used
1294 
1295    Level: intermediate
1296 
1297 .keywords: SNES, nonlinear, set, Gauss-Seidel
1298 
1299 .seealso: SNESSetGS(), SNESGetGS(), SNESSetUseGS()
1300 @*/
1301 
1302 PetscErrorCode SNESGetUseGS(SNES snes, PetscBool * usegs) {
1303   PetscFunctionBegin;
1304 
1305   *usegs = snes->usegs;
1306   PetscFunctionReturn(0);
1307 }
1308 
1309 #undef __FUNCT__
1310 #define __FUNCT__ "SNESSetGSSweeps"
1311 /*@
1312    SNESSetGSSweeps - Sets the number of sweeps of GS to use.
1313 
1314    Input Parameters:
1315 +  snes   - the SNES context
1316 -  sweeps  - the number of sweeps of GS to perform.
1317 
1318    Level: intermediate
1319 
1320 .keywords: SNES, nonlinear, set, Gauss-Siedel
1321 
1322 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGetGSSweeps()
1323 @*/
1324 
1325 PetscErrorCode SNESSetGSSweeps(SNES snes, PetscInt sweeps) {
1326   PetscFunctionBegin;
1327   snes->gssweeps = sweeps;
1328   PetscFunctionReturn(0);
1329 }
1330 
1331 
1332 #undef __FUNCT__
1333 #define __FUNCT__ "SNESGetGSSweeps"
1334 /*@
1335    SNESGetGSSweeps - Gets the number of sweeps GS will use.
1336 
1337    Input Parameters:
1338 .  snes   - the SNES context
1339 
1340    Output Parameters:
1341 .  sweeps  - the number of sweeps of GS to perform.
1342 
1343    Level: intermediate
1344 
1345 .keywords: SNES, nonlinear, set, Gauss-Siedel
1346 
1347 .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESSetGSSweeps()
1348 @*/
1349 PetscErrorCode SNESGetGSSweeps(SNES snes, PetscInt * sweeps) {
1350   PetscFunctionBegin;
1351   *sweeps = snes->gssweeps;
1352   PetscFunctionReturn(0);
1353 }
1354 
1355 #undef __FUNCT__
1356 #define __FUNCT__ "SNESPicardComputeFunction"
1357 PetscErrorCode  SNESPicardComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
1358 {
1359   PetscErrorCode ierr;
1360   PetscFunctionBegin;
1361   /*  A(x)*x - b(x) */
1362   ierr = (*snes->ops->computepjacobian)(snes,x,&snes->jacobian,&snes->jacobian_pre,&snes->matstruct,snes->jacP);CHKERRQ(ierr);
1363   ierr = (*snes->ops->computepfunction)(snes,x,f,snes->funP);CHKERRQ(ierr);
1364   ierr = VecView(x,PETSC_VIEWER_BINARY_WORLD);CHKERRQ(ierr);
1365   ierr = VecView(f,PETSC_VIEWER_BINARY_WORLD);CHKERRQ(ierr);
1366   ierr = VecScale(f,-1.0);CHKERRQ(ierr);
1367   ierr = MatMultAdd(snes->jacobian_pre,x,f,f);CHKERRQ(ierr);
1368   PetscFunctionReturn(0);
1369 }
1370 
1371 #undef __FUNCT__
1372 #define __FUNCT__ "SNESPicardComputeJacobian"
1373 PetscErrorCode  SNESPicardComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
1374 {
1375   PetscFunctionBegin;
1376   *flag = snes->matstruct;
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 #undef __FUNCT__
1381 #define __FUNCT__ "SNESSetPicard"
1382 /*@C
1383    SNESSetPicard - Use SNES to solve the semilinear-system A(x) x = b(x) via a Picard type iteration (Picard linearization)
1384 
1385    Logically Collective on SNES
1386 
1387    Input Parameters:
1388 +  snes - the SNES context
1389 .  r - vector to store function value
1390 .  func - function evaluation routine
1391 .  jmat - normally the same as mat but you can pass another matrix for which you compute the Jacobian of A(x) x - b(x) (see jmat below)
1392 .  mat - matrix to store A
1393 .  mfunc  - function to compute matrix value
1394 -  ctx - [optional] user-defined context for private data for the
1395          function evaluation routine (may be PETSC_NULL)
1396 
1397    Calling sequence of func:
1398 $    func (SNES snes,Vec x,Vec f,void *ctx);
1399 
1400 +  f - function vector
1401 -  ctx - optional user-defined function context
1402 
1403    Calling sequence of mfunc:
1404 $     mfunc (SNES snes,Vec x,Mat *jmat,Mat *mat,int *flag,void *ctx);
1405 
1406 +  x - input vector
1407 .  jmat - Form Jacobian matrix of A(x) x - b(x) if available, not there is really no reason to use it in this way since then you can just use SNESSetJacobian(),
1408           normally just pass mat in this location
1409 .  mat - form A(x) matrix
1410 .  flag - flag indicating information about the preconditioner matrix
1411    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1412 -  ctx - [optional] user-defined Jacobian context
1413 
1414    Notes:
1415     One can call SNESSetPicard() or SNESSetFunction() (and possibly SNESSetJacobian()) but cannot call both
1416 
1417 $     Solves the equation A(x) x = b(x) via the defect correction algorithm A(x^{n}) (x^{n+1} - x^{n}) = b(x^{n}) - A(x^{n})x^{n}
1418 $     Note that when an exact solver is used this corresponds to the "classic" Picard A(x^{n}) x^{n+1} = b(x^{n}) iteration.
1419 
1420      Run with -snes_mf_operator to solve the system with Newton's method using A(x^{n}) to construct the preconditioner.
1421 
1422    We implement the defect correction form of the Picard iteration because it converges much more generally when inexact linear solvers are used then
1423    the direct Picard iteration A(x^n) x^{n+1} = b(x^n)
1424 
1425    There is some controversity over the definition of a Picard iteration for nonlinear systems but almost everyone agrees that it involves a linear solve and some
1426    believe it is the iteration  A(x^{n}) x^{n+1} = b(x^{n}) hence we use the name Picard. If anyone has an authoritative  reference that defines the Picard iteration
1427    different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
1428 
1429    Level: beginner
1430 
1431 .keywords: SNES, nonlinear, set, function
1432 
1433 .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESGetPicard(), SNESLineSearchPreCheckPicard()
1434 @*/
1435 PetscErrorCode  SNESSetPicard(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),Mat jmat, Mat mat, PetscErrorCode (*mfunc)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1436 {
1437   PetscErrorCode ierr;
1438   PetscFunctionBegin;
1439   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1440   snes->ops->computepfunction = func;
1441   snes->ops->computepjacobian = mfunc;
1442   ierr = SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);CHKERRQ(ierr);
1443   ierr = SNESSetJacobian(snes,jmat,mat,SNESPicardComputeJacobian,ctx);CHKERRQ(ierr);
1444   PetscFunctionReturn(0);
1445 }
1446 
1447 #undef __FUNCT__
1448 #define __FUNCT__ "SNESSetComputeInitialGuess"
1449 /*@C
1450    SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
1451 
1452    Logically Collective on SNES
1453 
1454    Input Parameters:
1455 +  snes - the SNES context
1456 .  func - function evaluation routine
1457 -  ctx - [optional] user-defined context for private data for the
1458          function evaluation routine (may be PETSC_NULL)
1459 
1460    Calling sequence of func:
1461 $    func (SNES snes,Vec x,void *ctx);
1462 
1463 .  f - function vector
1464 -  ctx - optional user-defined function context
1465 
1466    Level: intermediate
1467 
1468 .keywords: SNES, nonlinear, set, function
1469 
1470 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
1471 @*/
1472 PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
1473 {
1474   PetscFunctionBegin;
1475   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1476   if (func) snes->ops->computeinitialguess = func;
1477   if (ctx)  snes->initialguessP            = ctx;
1478   PetscFunctionReturn(0);
1479 }
1480 
1481 /* --------------------------------------------------------------- */
1482 #undef __FUNCT__
1483 #define __FUNCT__ "SNESGetRhs"
1484 /*@C
1485    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
1486    it assumes a zero right hand side.
1487 
1488    Logically Collective on SNES
1489 
1490    Input Parameter:
1491 .  snes - the SNES context
1492 
1493    Output Parameter:
1494 .  rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null
1495 
1496    Level: intermediate
1497 
1498 .keywords: SNES, nonlinear, get, function, right hand side
1499 
1500 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
1501 @*/
1502 PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
1503 {
1504   PetscFunctionBegin;
1505   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1506   PetscValidPointer(rhs,2);
1507   *rhs = snes->vec_rhs;
1508   PetscFunctionReturn(0);
1509 }
1510 
1511 #undef __FUNCT__
1512 #define __FUNCT__ "SNESComputeFunction"
1513 /*@
1514    SNESComputeFunction - Calls the function that has been set with
1515                          SNESSetFunction().
1516 
1517    Collective on SNES
1518 
1519    Input Parameters:
1520 +  snes - the SNES context
1521 -  x - input vector
1522 
1523    Output Parameter:
1524 .  y - function vector, as set by SNESSetFunction()
1525 
1526    Notes:
1527    SNESComputeFunction() is typically used within nonlinear solvers
1528    implementations, so most users would not generally call this routine
1529    themselves.
1530 
1531    Level: developer
1532 
1533 .keywords: SNES, nonlinear, compute, function
1534 
1535 .seealso: SNESSetFunction(), SNESGetFunction()
1536 @*/
1537 PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
1538 {
1539   PetscErrorCode ierr;
1540 
1541   PetscFunctionBegin;
1542   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1543   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1544   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1545   PetscCheckSameComm(snes,1,x,2);
1546   PetscCheckSameComm(snes,1,y,3);
1547   VecValidValues(x,2,PETSC_TRUE);
1548 
1549   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1550   if (snes->ops->computefunction) {
1551     PetscStackPush("SNES user function");
1552     ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
1553     PetscStackPop;
1554   } else if (snes->dm) {
1555     ierr = DMComputeFunction(snes->dm,x,y);CHKERRQ(ierr);
1556   } else if (snes->vec_rhs) {
1557     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
1558   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
1559   if (snes->vec_rhs) {
1560     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
1561   }
1562   snes->nfuncs++;
1563   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1564   VecValidValues(y,3,PETSC_FALSE);
1565   PetscFunctionReturn(0);
1566 }
1567 
1568 #undef __FUNCT__
1569 #define __FUNCT__ "SNESComputeGS"
1570 /*@
1571    SNESComputeGS - Calls the Gauss-Seidel function that has been set with
1572                    SNESSetGS().
1573 
1574    Collective on SNES
1575 
1576    Input Parameters:
1577 +  snes - the SNES context
1578 .  x - input vector
1579 -  b - rhs vector
1580 
1581    Output Parameter:
1582 .  x - new solution vector
1583 
1584    Notes:
1585    SNESComputeGS() is typically used within composed nonlinear solver
1586    implementations, so most users would not generally call this routine
1587    themselves.
1588 
1589    Level: developer
1590 
1591 .keywords: SNES, nonlinear, compute, function
1592 
1593 .seealso: SNESSetGS(), SNESComputeFunction()
1594 @*/
1595 PetscErrorCode  SNESComputeGS(SNES snes,Vec b,Vec x)
1596 {
1597   PetscErrorCode ierr;
1598   PetscInt i;
1599 
1600   PetscFunctionBegin;
1601   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1602   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1603   if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,3);
1604   PetscCheckSameComm(snes,1,x,2);
1605   if(b) PetscCheckSameComm(snes,1,b,3);
1606   if (b) VecValidValues(b,2,PETSC_TRUE);
1607   ierr = PetscLogEventBegin(SNES_GSEval,snes,x,b,0);CHKERRQ(ierr);
1608   if (snes->ops->computegs) {
1609     for (i = 0; i < snes->gssweeps; i++) {
1610       PetscStackPush("SNES user GS");
1611       ierr = (*snes->ops->computegs)(snes,x,b,snes->gsP);CHKERRQ(ierr);
1612       PetscStackPop;
1613     }
1614   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetGS() before SNESComputeGS(), likely called from SNESSolve().");
1615   ierr = PetscLogEventEnd(SNES_GSEval,snes,x,b,0);CHKERRQ(ierr);
1616   VecValidValues(x,3,PETSC_FALSE);
1617   PetscFunctionReturn(0);
1618 }
1619 
1620 
1621 #undef __FUNCT__
1622 #define __FUNCT__ "SNESComputeJacobian"
1623 /*@
1624    SNESComputeJacobian - Computes the Jacobian matrix that has been
1625    set with SNESSetJacobian().
1626 
1627    Collective on SNES and Mat
1628 
1629    Input Parameters:
1630 +  snes - the SNES context
1631 -  x - input vector
1632 
1633    Output Parameters:
1634 +  A - Jacobian matrix
1635 .  B - optional preconditioning matrix
1636 -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1637 
1638   Options Database Keys:
1639 +    -snes_lag_preconditioner <lag>
1640 .    -snes_lag_jacobian <lag>
1641 .    -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
1642 .    -snes_compare_explicit_draw  - Compare the computed Jacobian to the finite difference Jacobian and draw the result
1643 .    -snes_compare_explicit_contour  - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
1644 .    -snes_compare_operator  - Make the comparison options above use the operator instead of the preconditioning matrix
1645 .    -snes_compare_coloring - Compute the finite differece Jacobian using coloring and display norms of difference
1646 .    -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
1647 .    -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
1648 .    -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
1649 .    -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
1650 .    -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
1651 -    -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences
1652 
1653 
1654    Notes:
1655    Most users should not need to explicitly call this routine, as it
1656    is used internally within the nonlinear solvers.
1657 
1658    See KSPSetOperators() for important information about setting the
1659    flag parameter.
1660 
1661    Level: developer
1662 
1663 .keywords: SNES, compute, Jacobian, matrix
1664 
1665 .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
1666 @*/
1667 PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1668 {
1669   PetscErrorCode ierr;
1670   PetscBool      flag;
1671 
1672   PetscFunctionBegin;
1673   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1674   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1675   PetscValidPointer(flg,5);
1676   PetscCheckSameComm(snes,1,X,2);
1677   VecValidValues(X,2,PETSC_TRUE);
1678   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1679 
1680   /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
1681 
1682   if (snes->lagjacobian == -2) {
1683     snes->lagjacobian = -1;
1684     ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr);
1685   } else if (snes->lagjacobian == -1) {
1686     *flg = SAME_PRECONDITIONER;
1687     ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr);
1688     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1689     if (flag) {
1690       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1691       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1692     }
1693     PetscFunctionReturn(0);
1694   } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) {
1695     *flg = SAME_PRECONDITIONER;
1696     ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr);
1697     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1698     if (flag) {
1699       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1700       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1701     }
1702     PetscFunctionReturn(0);
1703   }
1704 
1705   *flg = DIFFERENT_NONZERO_PATTERN;
1706   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1707   PetscStackPush("SNES user Jacobian function");
1708   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1709   PetscStackPop;
1710   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1711 
1712   if (snes->lagpreconditioner == -2) {
1713     ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr);
1714     snes->lagpreconditioner = -1;
1715   } else if (snes->lagpreconditioner == -1) {
1716     *flg = SAME_PRECONDITIONER;
1717     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1718   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1719     *flg = SAME_PRECONDITIONER;
1720     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1721   }
1722 
1723   /* make sure user returned a correct Jacobian and preconditioner */
1724   /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
1725     PetscValidHeaderSpecific(*B,MAT_CLASSID,4);   */
1726   {
1727     PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
1728     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit",&flag,PETSC_NULL);CHKERRQ(ierr);
1729     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",&flag_draw,PETSC_NULL);CHKERRQ(ierr);
1730     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",&flag_contour,PETSC_NULL);CHKERRQ(ierr);
1731     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_operator",&flag_operator,PETSC_NULL);CHKERRQ(ierr);
1732     if (flag || flag_draw || flag_contour) {
1733       Mat Bexp_mine = PETSC_NULL,Bexp,FDexp;
1734       MatStructure mstruct;
1735       PetscViewer vdraw,vstdout;
1736       PetscBool flg;
1737       if (flag_operator) {
1738         ierr = MatComputeExplicitOperator(*A,&Bexp_mine);CHKERRQ(ierr);
1739         Bexp = Bexp_mine;
1740       } else {
1741         /* See if the preconditioning matrix can be viewed and added directly */
1742         ierr = PetscTypeCompareAny((PetscObject)*B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr);
1743         if (flg) Bexp = *B;
1744         else {
1745           /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
1746           ierr = MatComputeExplicitOperator(*B,&Bexp_mine);CHKERRQ(ierr);
1747           Bexp = Bexp_mine;
1748         }
1749       }
1750       ierr = MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);CHKERRQ(ierr);
1751       ierr = SNESDefaultComputeJacobian(snes,X,&FDexp,&FDexp,&mstruct,NULL);CHKERRQ(ierr);
1752       ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&vstdout);CHKERRQ(ierr);
1753       if (flag_draw || flag_contour) {
1754         ierr = PetscViewerDrawOpen(((PetscObject)snes)->comm,0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr);
1755         if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);}
1756       } else vdraw = PETSC_NULL;
1757       ierr = PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator?"Jacobian":"preconditioning Jacobian");CHKERRQ(ierr);
1758       if (flag) {ierr = MatView(Bexp,vstdout);CHKERRQ(ierr);}
1759       if (vdraw) {ierr = MatView(Bexp,vdraw);CHKERRQ(ierr);}
1760       ierr = PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");CHKERRQ(ierr);
1761       if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);}
1762       if (vdraw) {ierr = MatView(FDexp,vdraw);CHKERRQ(ierr);}
1763       ierr = MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1764       ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");CHKERRQ(ierr);
1765       if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);}
1766       if (vdraw) {              /* Always use contour for the difference */
1767         ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);
1768         ierr = MatView(FDexp,vdraw);CHKERRQ(ierr);
1769         ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);
1770       }
1771       if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);}
1772       ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr);
1773       ierr = MatDestroy(&Bexp_mine);CHKERRQ(ierr);
1774       ierr = MatDestroy(&FDexp);CHKERRQ(ierr);
1775     }
1776   }
1777   {
1778     PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
1779     PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
1780     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring",&flag,PETSC_NULL);CHKERRQ(ierr);
1781     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_display",&flag_display,PETSC_NULL);CHKERRQ(ierr);
1782     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",&flag_draw,PETSC_NULL);CHKERRQ(ierr);
1783     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",&flag_contour,PETSC_NULL);CHKERRQ(ierr);
1784     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",&flag_threshold,PETSC_NULL);CHKERRQ(ierr);
1785     ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,PETSC_NULL);CHKERRQ(ierr);
1786     ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,PETSC_NULL);CHKERRQ(ierr);
1787     if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
1788       Mat Bfd;
1789       MatStructure mstruct;
1790       PetscViewer vdraw,vstdout;
1791       ISColoring iscoloring;
1792       MatFDColoring matfdcoloring;
1793       PetscErrorCode (*func)(SNES,Vec,Vec,void*);
1794       void *funcctx;
1795       PetscReal norm1,norm2,normmax;
1796 
1797       ierr = MatDuplicate(*B,MAT_DO_NOT_COPY_VALUES,&Bfd);CHKERRQ(ierr);
1798       ierr = MatGetColoring(Bfd,MATCOLORINGSL,&iscoloring);CHKERRQ(ierr);
1799       ierr = MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);CHKERRQ(ierr);
1800       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
1801 
1802       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
1803       ierr = SNESGetFunction(snes,PETSC_NULL,&func,&funcctx);CHKERRQ(ierr);
1804       ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode(*)(void))func,funcctx);CHKERRQ(ierr);
1805       ierr = PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);CHKERRQ(ierr);
1806       ierr = PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");CHKERRQ(ierr);
1807       ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr);
1808       ierr = MatFDColoringApply(Bfd,matfdcoloring,X,&mstruct,snes);CHKERRQ(ierr);
1809       ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr);
1810 
1811       ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&vstdout);CHKERRQ(ierr);
1812       if (flag_draw || flag_contour) {
1813         ierr = PetscViewerDrawOpen(((PetscObject)snes)->comm,0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr);
1814         if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);}
1815       } else vdraw = PETSC_NULL;
1816       ierr = PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");CHKERRQ(ierr);
1817       if (flag_display) {ierr = MatView(*B,vstdout);CHKERRQ(ierr);}
1818       if (vdraw) {ierr = MatView(*B,vdraw);CHKERRQ(ierr);}
1819       ierr = PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");CHKERRQ(ierr);
1820       if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);}
1821       if (vdraw) {ierr = MatView(Bfd,vdraw);CHKERRQ(ierr);}
1822       ierr = MatAYPX(Bfd,-1.0,*B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1823       ierr = MatNorm(Bfd,NORM_1,&norm1);CHKERRQ(ierr);
1824       ierr = MatNorm(Bfd,NORM_FROBENIUS,&norm2);CHKERRQ(ierr);
1825       ierr = MatNorm(Bfd,NORM_MAX,&normmax);CHKERRQ(ierr);
1826       ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%G normFrob=%G normmax=%G\n",norm1,norm2,normmax);CHKERRQ(ierr);
1827       if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);}
1828       if (vdraw) {              /* Always use contour for the difference */
1829         ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);
1830         ierr = MatView(Bfd,vdraw);CHKERRQ(ierr);
1831         ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);
1832       }
1833       if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);}
1834 
1835       if (flag_threshold) {
1836         PetscInt bs,rstart,rend,i;
1837         ierr = MatGetBlockSize(*B,&bs);CHKERRQ(ierr);
1838         ierr = MatGetOwnershipRange(*B,&rstart,&rend);CHKERRQ(ierr);
1839         for (i=rstart; i<rend; i++) {
1840           const PetscScalar *ba,*ca;
1841           const PetscInt *bj,*cj;
1842           PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
1843           PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
1844           ierr = MatGetRow(*B,i,&bn,&bj,&ba);CHKERRQ(ierr);
1845           ierr = MatGetRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr);
1846           if (bn != cn) SETERRQ(((PetscObject)*A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
1847           for (j=0; j<bn; j++) {
1848             PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
1849             if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
1850               maxentrycol = bj[j];
1851               maxentry = PetscRealPart(ba[j]);
1852             }
1853             if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
1854               maxdiffcol = bj[j];
1855               maxdiff = PetscRealPart(ca[j]);
1856             }
1857             if (rdiff > maxrdiff) {
1858               maxrdiffcol = bj[j];
1859               maxrdiff = rdiff;
1860             }
1861           }
1862           if (maxrdiff > 1) {
1863             ierr = PetscViewerASCIIPrintf(vstdout,"row %D (maxentry=%G at %D, maxdiff=%G at %D, maxrdiff=%G at %D):",i,maxentry,maxentrycol,maxdiff,maxdiffcol,maxrdiff,maxrdiffcol);CHKERRQ(ierr);
1864             for (j=0; j<bn; j++) {
1865               PetscReal rdiff;
1866               rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
1867               if (rdiff > 1) {
1868                 ierr = PetscViewerASCIIPrintf(vstdout," (%D,%G:%G)",bj[j],PetscRealPart(ba[j]),PetscRealPart(ca[j]));CHKERRQ(ierr);
1869               }
1870             }
1871             ierr = PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);CHKERRQ(ierr);
1872           }
1873           ierr = MatRestoreRow(*B,i,&bn,&bj,&ba);CHKERRQ(ierr);
1874           ierr = MatRestoreRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr);
1875         }
1876       }
1877       ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr);
1878       ierr = MatDestroy(&Bfd);CHKERRQ(ierr);
1879     }
1880   }
1881   PetscFunctionReturn(0);
1882 }
1883 
1884 #undef __FUNCT__
1885 #define __FUNCT__ "SNESSetJacobian"
1886 /*@C
1887    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1888    location to store the matrix.
1889 
1890    Logically Collective on SNES and Mat
1891 
1892    Input Parameters:
1893 +  snes - the SNES context
1894 .  A - Jacobian matrix
1895 .  B - preconditioner matrix (usually same as the Jacobian)
1896 .  func - Jacobian evaluation routine (if PETSC_NULL then SNES retains any previously set value)
1897 -  ctx - [optional] user-defined context for private data for the
1898          Jacobian evaluation routine (may be PETSC_NULL) (if PETSC_NULL then SNES retains any previously set value)
1899 
1900    Calling sequence of func:
1901 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1902 
1903 +  x - input vector
1904 .  A - Jacobian matrix
1905 .  B - preconditioner matrix, usually the same as A
1906 .  flag - flag indicating information about the preconditioner matrix
1907    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1908 -  ctx - [optional] user-defined Jacobian context
1909 
1910    Notes:
1911    See KSPSetOperators() for important information about setting the flag
1912    output parameter in the routine func().  Be sure to read this information!
1913 
1914    The routine func() takes Mat * as the matrix arguments rather than Mat.
1915    This allows the Jacobian evaluation routine to replace A and/or B with a
1916    completely new new matrix structure (not just different matrix elements)
1917    when appropriate, for instance, if the nonzero structure is changing
1918    throughout the global iterations.
1919 
1920    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
1921    each matrix.
1922 
1923    If using SNESDefaultComputeJacobianColor() to assemble a Jacobian, the ctx argument
1924    must be a MatFDColoring.
1925 
1926    Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian.  One common
1927    example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
1928 
1929    Level: beginner
1930 
1931 .keywords: SNES, nonlinear, set, Jacobian, matrix
1932 
1933 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
1934 @*/
1935 PetscErrorCode  SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1936 {
1937   PetscErrorCode ierr;
1938 
1939   PetscFunctionBegin;
1940   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1941   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
1942   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
1943   if (A) PetscCheckSameComm(snes,1,A,2);
1944   if (B) PetscCheckSameComm(snes,1,B,3);
1945   if (func) snes->ops->computejacobian = func;
1946   if (ctx)  snes->jacP                 = ctx;
1947   if (A) {
1948     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1949     ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr);
1950     snes->jacobian = A;
1951   }
1952   if (B) {
1953     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1954     ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr);
1955     snes->jacobian_pre = B;
1956   }
1957   PetscFunctionReturn(0);
1958 }
1959 
1960 #undef __FUNCT__
1961 #define __FUNCT__ "SNESGetJacobian"
1962 /*@C
1963    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1964    provided context for evaluating the Jacobian.
1965 
1966    Not Collective, but Mat object will be parallel if SNES object is
1967 
1968    Input Parameter:
1969 .  snes - the nonlinear solver context
1970 
1971    Output Parameters:
1972 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1973 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1974 .  func - location to put Jacobian function (or PETSC_NULL)
1975 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1976 
1977    Level: advanced
1978 
1979 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1980 @*/
1981 PetscErrorCode  SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1982 {
1983   PetscFunctionBegin;
1984   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1985   if (A)    *A    = snes->jacobian;
1986   if (B)    *B    = snes->jacobian_pre;
1987   if (func) *func = snes->ops->computejacobian;
1988   if (ctx)  *ctx  = snes->jacP;
1989   PetscFunctionReturn(0);
1990 }
1991 
1992 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1993 
1994 #undef __FUNCT__
1995 #define __FUNCT__ "SNESSetUp"
1996 /*@
1997    SNESSetUp - Sets up the internal data structures for the later use
1998    of a nonlinear solver.
1999 
2000    Collective on SNES
2001 
2002    Input Parameters:
2003 .  snes - the SNES context
2004 
2005    Notes:
2006    For basic use of the SNES solvers the user need not explicitly call
2007    SNESSetUp(), since these actions will automatically occur during
2008    the call to SNESSolve().  However, if one wishes to control this
2009    phase separately, SNESSetUp() should be called after SNESCreate()
2010    and optional routines of the form SNESSetXXX(), but before SNESSolve().
2011 
2012    Level: advanced
2013 
2014 .keywords: SNES, nonlinear, setup
2015 
2016 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
2017 @*/
2018 PetscErrorCode  SNESSetUp(SNES snes)
2019 {
2020   PetscErrorCode ierr;
2021 
2022   PetscFunctionBegin;
2023   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2024   if (snes->setupcalled) PetscFunctionReturn(0);
2025 
2026   if (!((PetscObject)snes)->type_name) {
2027     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
2028   }
2029 
2030   ierr = SNESGetFunction(snes,&snes->vec_func,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2031   if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
2032   if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
2033 
2034   if (!snes->vec_sol_update /* && snes->vec_sol */) {
2035     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
2036     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
2037   }
2038 
2039   if (!snes->ops->computejacobian && snes->dm) {
2040     Mat J,B;
2041     ierr = DMCreateMatrix(snes->dm,MATAIJ,&B);CHKERRQ(ierr);
2042     if (snes->mf_operator) {
2043       ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
2044       ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr);
2045       ierr = MatSetFromOptions(J);CHKERRQ(ierr);
2046     } else {
2047       J = B;
2048       ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr);
2049     }
2050     ierr = SNESSetJacobian(snes,J,B,SNESDMComputeJacobian,PETSC_NULL);CHKERRQ(ierr);
2051     ierr = MatDestroy(&J);CHKERRQ(ierr);
2052     ierr = MatDestroy(&B);CHKERRQ(ierr);
2053   } else if (!snes->jacobian && snes->ops->computejacobian == MatMFFDComputeJacobian) {
2054     Mat J;
2055     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
2056     ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr);
2057     ierr = MatSetFromOptions(J);CHKERRQ(ierr);
2058     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
2059     ierr = MatDestroy(&J);CHKERRQ(ierr);
2060   } else if (snes->dm && snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
2061     Mat J,B;
2062     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
2063     ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr);
2064     ierr = MatSetFromOptions(J);CHKERRQ(ierr);
2065     ierr = DMCreateMatrix(snes->dm,MATAIJ,&B);CHKERRQ(ierr);
2066     ierr = SNESSetJacobian(snes,J,B,SNESDMComputeJacobian,snes->funP);CHKERRQ(ierr);
2067     ierr = MatDestroy(&J);CHKERRQ(ierr);
2068     ierr = MatDestroy(&B);CHKERRQ(ierr);
2069   } else if (snes->dm && !snes->jacobian_pre){
2070     Mat J;
2071     ierr = DMCreateMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr);
2072     ierr = SNESSetJacobian(snes,J,J,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2073     ierr = MatDestroy(&J);CHKERRQ(ierr);
2074   }
2075   if (!snes->ops->computefunction && !snes->dm) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must provide a residual function with SNESSetFunction() or call SNESSetDM()");
2076   if (!snes->vec_func) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must provide a vector when calling SNESSetFunction() or call SNESSetDM()");
2077 
2078   if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);}
2079 
2080   if (snes->ops->usercompute && !snes->user) {
2081     ierr = (*snes->ops->usercompute)(snes,(void**)&snes->user);CHKERRQ(ierr);
2082   }
2083 
2084   if (snes->ops->setup) {
2085     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
2086   }
2087 
2088   snes->setupcalled = PETSC_TRUE;
2089   PetscFunctionReturn(0);
2090 }
2091 
2092 #undef __FUNCT__
2093 #define __FUNCT__ "SNESReset"
2094 /*@
2095    SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
2096 
2097    Collective on SNES
2098 
2099    Input Parameter:
2100 .  snes - iterative context obtained from SNESCreate()
2101 
2102    Level: intermediate
2103 
2104    Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
2105 
2106 .keywords: SNES, destroy
2107 
2108 .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2109 @*/
2110 PetscErrorCode  SNESReset(SNES snes)
2111 {
2112   PetscErrorCode ierr;
2113 
2114   PetscFunctionBegin;
2115   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2116   if (snes->ops->userdestroy && snes->user) {
2117     ierr       = (*snes->ops->userdestroy)((void**)&snes->user);CHKERRQ(ierr);
2118     snes->user = PETSC_NULL;
2119   }
2120   if (snes->pc) {
2121     ierr = SNESReset(snes->pc);CHKERRQ(ierr);
2122   }
2123 
2124   if (snes->ops->reset) {
2125     ierr = (*snes->ops->reset)(snes);CHKERRQ(ierr);
2126   }
2127   if (snes->ksp) {ierr = KSPReset(snes->ksp);CHKERRQ(ierr);}
2128   ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr);
2129   ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr);
2130   ierr = VecDestroy(&snes->vec_sol_update);CHKERRQ(ierr);
2131   ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr);
2132   ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr);
2133   ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr);
2134   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2135   if (snes->vwork) {ierr = VecDestroyVecs(snes->nvwork,&snes->vwork);CHKERRQ(ierr);}
2136   snes->nwork = snes->nvwork = 0;
2137   snes->setupcalled = PETSC_FALSE;
2138   PetscFunctionReturn(0);
2139 }
2140 
2141 #undef __FUNCT__
2142 #define __FUNCT__ "SNESDestroy"
2143 /*@
2144    SNESDestroy - Destroys the nonlinear solver context that was created
2145    with SNESCreate().
2146 
2147    Collective on SNES
2148 
2149    Input Parameter:
2150 .  snes - the SNES context
2151 
2152    Level: beginner
2153 
2154 .keywords: SNES, nonlinear, destroy
2155 
2156 .seealso: SNESCreate(), SNESSolve()
2157 @*/
2158 PetscErrorCode  SNESDestroy(SNES *snes)
2159 {
2160   PetscErrorCode ierr;
2161 
2162   PetscFunctionBegin;
2163   if (!*snes) PetscFunctionReturn(0);
2164   PetscValidHeaderSpecific((*snes),SNES_CLASSID,1);
2165   if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; PetscFunctionReturn(0);}
2166 
2167   ierr = SNESReset((*snes));CHKERRQ(ierr);
2168   ierr = SNESDestroy(&(*snes)->pc);CHKERRQ(ierr);
2169 
2170   /* if memory was published with AMS then destroy it */
2171   ierr = PetscObjectDepublish((*snes));CHKERRQ(ierr);
2172   if ((*snes)->ops->destroy) {ierr = (*((*snes))->ops->destroy)((*snes));CHKERRQ(ierr);}
2173 
2174   ierr = DMDestroy(&(*snes)->dm);CHKERRQ(ierr);
2175   ierr = KSPDestroy(&(*snes)->ksp);CHKERRQ(ierr);
2176 
2177   ierr = PetscFree((*snes)->kspconvctx);CHKERRQ(ierr);
2178   if ((*snes)->ops->convergeddestroy) {
2179     ierr = (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);CHKERRQ(ierr);
2180   }
2181   if ((*snes)->conv_malloc) {
2182     ierr = PetscFree((*snes)->conv_hist);CHKERRQ(ierr);
2183     ierr = PetscFree((*snes)->conv_hist_its);CHKERRQ(ierr);
2184   }
2185   ierr = PetscViewerDestroy(&(*snes)->ls_monitor);CHKERRQ(ierr);
2186   ierr = SNESMonitorCancel((*snes));CHKERRQ(ierr);
2187   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
2188  PetscFunctionReturn(0);
2189 }
2190 
2191 /* ----------- Routines to set solver parameters ---------- */
2192 
2193 #undef __FUNCT__
2194 #define __FUNCT__ "SNESSetLagPreconditioner"
2195 /*@
2196    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
2197 
2198    Logically Collective on SNES
2199 
2200    Input Parameters:
2201 +  snes - the SNES context
2202 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2203          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2204 
2205    Options Database Keys:
2206 .    -snes_lag_preconditioner <lag>
2207 
2208    Notes:
2209    The default is 1
2210    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2211    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
2212 
2213    Level: intermediate
2214 
2215 .keywords: SNES, nonlinear, set, convergence, tolerances
2216 
2217 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2218 
2219 @*/
2220 PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2221 {
2222   PetscFunctionBegin;
2223   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2224   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2225   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2226   PetscValidLogicalCollectiveInt(snes,lag,2);
2227   snes->lagpreconditioner = lag;
2228   PetscFunctionReturn(0);
2229 }
2230 
2231 #undef __FUNCT__
2232 #define __FUNCT__ "SNESSetGridSequence"
2233 /*@
2234    SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
2235 
2236    Logically Collective on SNES
2237 
2238    Input Parameters:
2239 +  snes - the SNES context
2240 -  steps - the number of refinements to do, defaults to 0
2241 
2242    Options Database Keys:
2243 .    -snes_grid_sequence <steps>
2244 
2245    Level: intermediate
2246 
2247    Notes:
2248    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
2249 
2250 .keywords: SNES, nonlinear, set, convergence, tolerances
2251 
2252 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2253 
2254 @*/
2255 PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
2256 {
2257   PetscFunctionBegin;
2258   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2259   PetscValidLogicalCollectiveInt(snes,steps,2);
2260   snes->gridsequence = steps;
2261   PetscFunctionReturn(0);
2262 }
2263 
2264 #undef __FUNCT__
2265 #define __FUNCT__ "SNESGetLagPreconditioner"
2266 /*@
2267    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
2268 
2269    Not Collective
2270 
2271    Input Parameter:
2272 .  snes - the SNES context
2273 
2274    Output Parameter:
2275 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2276          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2277 
2278    Options Database Keys:
2279 .    -snes_lag_preconditioner <lag>
2280 
2281    Notes:
2282    The default is 1
2283    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2284 
2285    Level: intermediate
2286 
2287 .keywords: SNES, nonlinear, set, convergence, tolerances
2288 
2289 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
2290 
2291 @*/
2292 PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2293 {
2294   PetscFunctionBegin;
2295   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2296   *lag = snes->lagpreconditioner;
2297   PetscFunctionReturn(0);
2298 }
2299 
2300 #undef __FUNCT__
2301 #define __FUNCT__ "SNESSetLagJacobian"
2302 /*@
2303    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
2304      often the preconditioner is rebuilt.
2305 
2306    Logically Collective on SNES
2307 
2308    Input Parameters:
2309 +  snes - the SNES context
2310 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2311          the Jacobian is built etc. -2 means rebuild at next chance but then never again
2312 
2313    Options Database Keys:
2314 .    -snes_lag_jacobian <lag>
2315 
2316    Notes:
2317    The default is 1
2318    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2319    If  -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used, use -2 to indicate you want it recomputed
2320    at the next Newton step but never again (unless it is reset to another value)
2321 
2322    Level: intermediate
2323 
2324 .keywords: SNES, nonlinear, set, convergence, tolerances
2325 
2326 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
2327 
2328 @*/
2329 PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
2330 {
2331   PetscFunctionBegin;
2332   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2333   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2334   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2335   PetscValidLogicalCollectiveInt(snes,lag,2);
2336   snes->lagjacobian = lag;
2337   PetscFunctionReturn(0);
2338 }
2339 
2340 #undef __FUNCT__
2341 #define __FUNCT__ "SNESGetLagJacobian"
2342 /*@
2343    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
2344 
2345    Not Collective
2346 
2347    Input Parameter:
2348 .  snes - the SNES context
2349 
2350    Output Parameter:
2351 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2352          the Jacobian is built etc.
2353 
2354    Options Database Keys:
2355 .    -snes_lag_jacobian <lag>
2356 
2357    Notes:
2358    The default is 1
2359    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2360 
2361    Level: intermediate
2362 
2363 .keywords: SNES, nonlinear, set, convergence, tolerances
2364 
2365 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
2366 
2367 @*/
2368 PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
2369 {
2370   PetscFunctionBegin;
2371   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2372   *lag = snes->lagjacobian;
2373   PetscFunctionReturn(0);
2374 }
2375 
2376 #undef __FUNCT__
2377 #define __FUNCT__ "SNESSetTolerances"
2378 /*@
2379    SNESSetTolerances - Sets various parameters used in convergence tests.
2380 
2381    Logically Collective on SNES
2382 
2383    Input Parameters:
2384 +  snes - the SNES context
2385 .  abstol - absolute convergence tolerance
2386 .  rtol - relative convergence tolerance
2387 .  stol -  convergence tolerance in terms of the norm
2388            of the change in the solution between steps
2389 .  maxit - maximum number of iterations
2390 -  maxf - maximum number of function evaluations
2391 
2392    Options Database Keys:
2393 +    -snes_atol <abstol> - Sets abstol
2394 .    -snes_rtol <rtol> - Sets rtol
2395 .    -snes_stol <stol> - Sets stol
2396 .    -snes_max_it <maxit> - Sets maxit
2397 -    -snes_max_funcs <maxf> - Sets maxf
2398 
2399    Notes:
2400    The default maximum number of iterations is 50.
2401    The default maximum number of function evaluations is 1000.
2402 
2403    Level: intermediate
2404 
2405 .keywords: SNES, nonlinear, set, convergence, tolerances
2406 
2407 .seealso: SNESSetTrustRegionTolerance()
2408 @*/
2409 PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
2410 {
2411   PetscFunctionBegin;
2412   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2413   PetscValidLogicalCollectiveReal(snes,abstol,2);
2414   PetscValidLogicalCollectiveReal(snes,rtol,3);
2415   PetscValidLogicalCollectiveReal(snes,stol,4);
2416   PetscValidLogicalCollectiveInt(snes,maxit,5);
2417   PetscValidLogicalCollectiveInt(snes,maxf,6);
2418 
2419   if (abstol != PETSC_DEFAULT) {
2420     if (abstol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol);
2421     snes->abstol = abstol;
2422   }
2423   if (rtol != PETSC_DEFAULT) {
2424     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %G must be non-negative and less than 1.0",rtol);
2425     snes->rtol = rtol;
2426   }
2427   if (stol != PETSC_DEFAULT) {
2428     if (stol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %G must be non-negative",stol);
2429     snes->xtol = stol;
2430   }
2431   if (maxit != PETSC_DEFAULT) {
2432     if (maxit < 0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
2433     snes->max_its = maxit;
2434   }
2435   if (maxf != PETSC_DEFAULT) {
2436     if (maxf < 0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
2437     snes->max_funcs = maxf;
2438   }
2439   PetscFunctionReturn(0);
2440 }
2441 
2442 #undef __FUNCT__
2443 #define __FUNCT__ "SNESGetTolerances"
2444 /*@
2445    SNESGetTolerances - Gets various parameters used in convergence tests.
2446 
2447    Not Collective
2448 
2449    Input Parameters:
2450 +  snes - the SNES context
2451 .  atol - absolute convergence tolerance
2452 .  rtol - relative convergence tolerance
2453 .  stol -  convergence tolerance in terms of the norm
2454            of the change in the solution between steps
2455 .  maxit - maximum number of iterations
2456 -  maxf - maximum number of function evaluations
2457 
2458    Notes:
2459    The user can specify PETSC_NULL for any parameter that is not needed.
2460 
2461    Level: intermediate
2462 
2463 .keywords: SNES, nonlinear, get, convergence, tolerances
2464 
2465 .seealso: SNESSetTolerances()
2466 @*/
2467 PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
2468 {
2469   PetscFunctionBegin;
2470   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2471   if (atol)  *atol  = snes->abstol;
2472   if (rtol)  *rtol  = snes->rtol;
2473   if (stol)  *stol  = snes->xtol;
2474   if (maxit) *maxit = snes->max_its;
2475   if (maxf)  *maxf  = snes->max_funcs;
2476   PetscFunctionReturn(0);
2477 }
2478 
2479 #undef __FUNCT__
2480 #define __FUNCT__ "SNESSetTrustRegionTolerance"
2481 /*@
2482    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
2483 
2484    Logically Collective on SNES
2485 
2486    Input Parameters:
2487 +  snes - the SNES context
2488 -  tol - tolerance
2489 
2490    Options Database Key:
2491 .  -snes_trtol <tol> - Sets tol
2492 
2493    Level: intermediate
2494 
2495 .keywords: SNES, nonlinear, set, trust region, tolerance
2496 
2497 .seealso: SNESSetTolerances()
2498 @*/
2499 PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
2500 {
2501   PetscFunctionBegin;
2502   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2503   PetscValidLogicalCollectiveReal(snes,tol,2);
2504   snes->deltatol = tol;
2505   PetscFunctionReturn(0);
2506 }
2507 
2508 /*
2509    Duplicate the lg monitors for SNES from KSP; for some reason with
2510    dynamic libraries things don't work under Sun4 if we just use
2511    macros instead of functions
2512 */
2513 #undef __FUNCT__
2514 #define __FUNCT__ "SNESMonitorLG"
2515 PetscErrorCode  SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
2516 {
2517   PetscErrorCode ierr;
2518 
2519   PetscFunctionBegin;
2520   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2521   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
2522   PetscFunctionReturn(0);
2523 }
2524 
2525 #undef __FUNCT__
2526 #define __FUNCT__ "SNESMonitorLGCreate"
2527 PetscErrorCode  SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2528 {
2529   PetscErrorCode ierr;
2530 
2531   PetscFunctionBegin;
2532   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
2533   PetscFunctionReturn(0);
2534 }
2535 
2536 #undef __FUNCT__
2537 #define __FUNCT__ "SNESMonitorLGDestroy"
2538 PetscErrorCode  SNESMonitorLGDestroy(PetscDrawLG *draw)
2539 {
2540   PetscErrorCode ierr;
2541 
2542   PetscFunctionBegin;
2543   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
2544   PetscFunctionReturn(0);
2545 }
2546 
2547 extern PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
2548 #undef __FUNCT__
2549 #define __FUNCT__ "SNESMonitorLGRange"
2550 PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
2551 {
2552   PetscDrawLG      lg;
2553   PetscErrorCode   ierr;
2554   PetscReal        x,y,per;
2555   PetscViewer      v = (PetscViewer)monctx;
2556   static PetscReal prev; /* should be in the context */
2557   PetscDraw        draw;
2558   PetscFunctionBegin;
2559   if (!monctx) {
2560     MPI_Comm    comm;
2561 
2562     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
2563     v      = PETSC_VIEWER_DRAW_(comm);
2564   }
2565   ierr   = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
2566   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2567   ierr   = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2568   ierr   = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
2569   x = (PetscReal) n;
2570   if (rnorm > 0.0) y = log10(rnorm); else y = -15.0;
2571   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2572   if (n < 20 || !(n % 5)) {
2573     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2574   }
2575 
2576   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
2577   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2578   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2579   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
2580   ierr =  SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr);
2581   x = (PetscReal) n;
2582   y = 100.0*per;
2583   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2584   if (n < 20 || !(n % 5)) {
2585     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2586   }
2587 
2588   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
2589   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2590   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2591   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
2592   x = (PetscReal) n;
2593   y = (prev - rnorm)/prev;
2594   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2595   if (n < 20 || !(n % 5)) {
2596     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2597   }
2598 
2599   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
2600   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2601   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2602   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
2603   x = (PetscReal) n;
2604   y = (prev - rnorm)/(prev*per);
2605   if (n > 2) { /*skip initial crazy value */
2606     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2607   }
2608   if (n < 20 || !(n % 5)) {
2609     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2610   }
2611   prev = rnorm;
2612   PetscFunctionReturn(0);
2613 }
2614 
2615 #undef __FUNCT__
2616 #define __FUNCT__ "SNESMonitorLGRangeCreate"
2617 PetscErrorCode  SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2618 {
2619   PetscErrorCode ierr;
2620 
2621   PetscFunctionBegin;
2622   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
2623   PetscFunctionReturn(0);
2624 }
2625 
2626 #undef __FUNCT__
2627 #define __FUNCT__ "SNESMonitorLGRangeDestroy"
2628 PetscErrorCode  SNESMonitorLGRangeDestroy(PetscDrawLG *draw)
2629 {
2630   PetscErrorCode ierr;
2631 
2632   PetscFunctionBegin;
2633   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
2634   PetscFunctionReturn(0);
2635 }
2636 
2637 #undef __FUNCT__
2638 #define __FUNCT__ "SNESMonitor"
2639 /*@
2640    SNESMonitor - runs the user provided monitor routines, if they exist
2641 
2642    Collective on SNES
2643 
2644    Input Parameters:
2645 +  snes - nonlinear solver context obtained from SNESCreate()
2646 .  iter - iteration number
2647 -  rnorm - relative norm of the residual
2648 
2649    Notes:
2650    This routine is called by the SNES implementations.
2651    It does not typically need to be called by the user.
2652 
2653    Level: developer
2654 
2655 .seealso: SNESMonitorSet()
2656 @*/
2657 PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
2658 {
2659   PetscErrorCode ierr;
2660   PetscInt       i,n = snes->numbermonitors;
2661 
2662   PetscFunctionBegin;
2663   for (i=0; i<n; i++) {
2664     ierr = (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);CHKERRQ(ierr);
2665   }
2666   PetscFunctionReturn(0);
2667 }
2668 
2669 /* ------------ Routines to set performance monitoring options ----------- */
2670 
2671 #undef __FUNCT__
2672 #define __FUNCT__ "SNESMonitorSet"
2673 /*@C
2674    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
2675    iteration of the nonlinear solver to display the iteration's
2676    progress.
2677 
2678    Logically Collective on SNES
2679 
2680    Input Parameters:
2681 +  snes - the SNES context
2682 .  func - monitoring routine
2683 .  mctx - [optional] user-defined context for private data for the
2684           monitor routine (use PETSC_NULL if no context is desired)
2685 -  monitordestroy - [optional] routine that frees monitor context
2686           (may be PETSC_NULL)
2687 
2688    Calling sequence of func:
2689 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
2690 
2691 +    snes - the SNES context
2692 .    its - iteration number
2693 .    norm - 2-norm function value (may be estimated)
2694 -    mctx - [optional] monitoring context
2695 
2696    Options Database Keys:
2697 +    -snes_monitor        - sets SNESMonitorDefault()
2698 .    -snes_monitor_draw    - sets line graph monitor,
2699                             uses SNESMonitorLGCreate()
2700 -    -snes_monitor_cancel - cancels all monitors that have
2701                             been hardwired into a code by
2702                             calls to SNESMonitorSet(), but
2703                             does not cancel those set via
2704                             the options database.
2705 
2706    Notes:
2707    Several different monitoring routines may be set by calling
2708    SNESMonitorSet() multiple times; all will be called in the
2709    order in which they were set.
2710 
2711    Fortran notes: Only a single monitor function can be set for each SNES object
2712 
2713    Level: intermediate
2714 
2715 .keywords: SNES, nonlinear, set, monitor
2716 
2717 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
2718 @*/
2719 PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
2720 {
2721   PetscInt       i;
2722   PetscErrorCode ierr;
2723 
2724   PetscFunctionBegin;
2725   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2726   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
2727   for (i=0; i<snes->numbermonitors;i++) {
2728     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
2729       if (monitordestroy) {
2730         ierr = (*monitordestroy)(&mctx);CHKERRQ(ierr);
2731       }
2732       PetscFunctionReturn(0);
2733     }
2734   }
2735   snes->monitor[snes->numbermonitors]           = monitor;
2736   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
2737   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
2738   PetscFunctionReturn(0);
2739 }
2740 
2741 #undef __FUNCT__
2742 #define __FUNCT__ "SNESMonitorCancel"
2743 /*@C
2744    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
2745 
2746    Logically Collective on SNES
2747 
2748    Input Parameters:
2749 .  snes - the SNES context
2750 
2751    Options Database Key:
2752 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
2753     into a code by calls to SNESMonitorSet(), but does not cancel those
2754     set via the options database
2755 
2756    Notes:
2757    There is no way to clear one specific monitor from a SNES object.
2758 
2759    Level: intermediate
2760 
2761 .keywords: SNES, nonlinear, set, monitor
2762 
2763 .seealso: SNESMonitorDefault(), SNESMonitorSet()
2764 @*/
2765 PetscErrorCode  SNESMonitorCancel(SNES snes)
2766 {
2767   PetscErrorCode ierr;
2768   PetscInt       i;
2769 
2770   PetscFunctionBegin;
2771   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2772   for (i=0; i<snes->numbermonitors; i++) {
2773     if (snes->monitordestroy[i]) {
2774       ierr = (*snes->monitordestroy[i])(&snes->monitorcontext[i]);CHKERRQ(ierr);
2775     }
2776   }
2777   snes->numbermonitors = 0;
2778   PetscFunctionReturn(0);
2779 }
2780 
2781 #undef __FUNCT__
2782 #define __FUNCT__ "SNESSetConvergenceTest"
2783 /*@C
2784    SNESSetConvergenceTest - Sets the function that is to be used
2785    to test for convergence of the nonlinear iterative solution.
2786 
2787    Logically Collective on SNES
2788 
2789    Input Parameters:
2790 +  snes - the SNES context
2791 .  func - routine to test for convergence
2792 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
2793 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
2794 
2795    Calling sequence of func:
2796 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
2797 
2798 +    snes - the SNES context
2799 .    it - current iteration (0 is the first and is before any Newton step)
2800 .    cctx - [optional] convergence context
2801 .    reason - reason for convergence/divergence
2802 .    xnorm - 2-norm of current iterate
2803 .    gnorm - 2-norm of current step
2804 -    f - 2-norm of function
2805 
2806    Level: advanced
2807 
2808 .keywords: SNES, nonlinear, set, convergence, test
2809 
2810 .seealso: SNESDefaultConverged(), SNESSkipConverged()
2811 @*/
2812 PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2813 {
2814   PetscErrorCode ierr;
2815 
2816   PetscFunctionBegin;
2817   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2818   if (!func) func = SNESSkipConverged;
2819   if (snes->ops->convergeddestroy) {
2820     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
2821   }
2822   snes->ops->converged        = func;
2823   snes->ops->convergeddestroy = destroy;
2824   snes->cnvP                  = cctx;
2825   PetscFunctionReturn(0);
2826 }
2827 
2828 #undef __FUNCT__
2829 #define __FUNCT__ "SNESGetConvergedReason"
2830 /*@
2831    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
2832 
2833    Not Collective
2834 
2835    Input Parameter:
2836 .  snes - the SNES context
2837 
2838    Output Parameter:
2839 .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
2840             manual pages for the individual convergence tests for complete lists
2841 
2842    Level: intermediate
2843 
2844    Notes: Can only be called after the call the SNESSolve() is complete.
2845 
2846 .keywords: SNES, nonlinear, set, convergence, test
2847 
2848 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
2849 @*/
2850 PetscErrorCode  SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
2851 {
2852   PetscFunctionBegin;
2853   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2854   PetscValidPointer(reason,2);
2855   *reason = snes->reason;
2856   PetscFunctionReturn(0);
2857 }
2858 
2859 #undef __FUNCT__
2860 #define __FUNCT__ "SNESSetConvergenceHistory"
2861 /*@
2862    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
2863 
2864    Logically Collective on SNES
2865 
2866    Input Parameters:
2867 +  snes - iterative context obtained from SNESCreate()
2868 .  a   - array to hold history, this array will contain the function norms computed at each step
2869 .  its - integer array holds the number of linear iterations for each solve.
2870 .  na  - size of a and its
2871 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
2872            else it continues storing new values for new nonlinear solves after the old ones
2873 
2874    Notes:
2875    If 'a' and 'its' are PETSC_NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
2876    default array of length 10000 is allocated.
2877 
2878    This routine is useful, e.g., when running a code for purposes
2879    of accurate performance monitoring, when no I/O should be done
2880    during the section of code that is being timed.
2881 
2882    Level: intermediate
2883 
2884 .keywords: SNES, set, convergence, history
2885 
2886 .seealso: SNESGetConvergenceHistory()
2887 
2888 @*/
2889 PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool  reset)
2890 {
2891   PetscErrorCode ierr;
2892 
2893   PetscFunctionBegin;
2894   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2895   if (na)  PetscValidScalarPointer(a,2);
2896   if (its) PetscValidIntPointer(its,3);
2897   if (na == PETSC_DECIDE || na == PETSC_DEFAULT || !a) {
2898     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2899     ierr = PetscMalloc(na*sizeof(PetscReal),&a);CHKERRQ(ierr);
2900     ierr = PetscMalloc(na*sizeof(PetscInt),&its);CHKERRQ(ierr);
2901     snes->conv_malloc   = PETSC_TRUE;
2902   }
2903   snes->conv_hist       = a;
2904   snes->conv_hist_its   = its;
2905   snes->conv_hist_max   = na;
2906   snes->conv_hist_len   = 0;
2907   snes->conv_hist_reset = reset;
2908   PetscFunctionReturn(0);
2909 }
2910 
2911 #if defined(PETSC_HAVE_MATLAB_ENGINE)
2912 #include <engine.h>   /* MATLAB include file */
2913 #include <mex.h>      /* MATLAB include file */
2914 EXTERN_C_BEGIN
2915 #undef __FUNCT__
2916 #define __FUNCT__ "SNESGetConvergenceHistoryMatlab"
2917 mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
2918 {
2919   mxArray        *mat;
2920   PetscInt       i;
2921   PetscReal      *ar;
2922 
2923   PetscFunctionBegin;
2924   mat  = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
2925   ar   = (PetscReal*) mxGetData(mat);
2926   for (i=0; i<snes->conv_hist_len; i++) {
2927     ar[i] = snes->conv_hist[i];
2928   }
2929   PetscFunctionReturn(mat);
2930 }
2931 EXTERN_C_END
2932 #endif
2933 
2934 
2935 #undef __FUNCT__
2936 #define __FUNCT__ "SNESGetConvergenceHistory"
2937 /*@C
2938    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
2939 
2940    Not Collective
2941 
2942    Input Parameter:
2943 .  snes - iterative context obtained from SNESCreate()
2944 
2945    Output Parameters:
2946 .  a   - array to hold history
2947 .  its - integer array holds the number of linear iterations (or
2948          negative if not converged) for each solve.
2949 -  na  - size of a and its
2950 
2951    Notes:
2952     The calling sequence for this routine in Fortran is
2953 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
2954 
2955    This routine is useful, e.g., when running a code for purposes
2956    of accurate performance monitoring, when no I/O should be done
2957    during the section of code that is being timed.
2958 
2959    Level: intermediate
2960 
2961 .keywords: SNES, get, convergence, history
2962 
2963 .seealso: SNESSetConvergencHistory()
2964 
2965 @*/
2966 PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
2967 {
2968   PetscFunctionBegin;
2969   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2970   if (a)   *a   = snes->conv_hist;
2971   if (its) *its = snes->conv_hist_its;
2972   if (na)  *na  = snes->conv_hist_len;
2973   PetscFunctionReturn(0);
2974 }
2975 
2976 #undef __FUNCT__
2977 #define __FUNCT__ "SNESSetUpdate"
2978 /*@C
2979   SNESSetUpdate - Sets the general-purpose update function called
2980   at the beginning of every iteration of the nonlinear solve. Specifically
2981   it is called just before the Jacobian is "evaluated".
2982 
2983   Logically Collective on SNES
2984 
2985   Input Parameters:
2986 . snes - The nonlinear solver context
2987 . func - The function
2988 
2989   Calling sequence of func:
2990 . func (SNES snes, PetscInt step);
2991 
2992 . step - The current step of the iteration
2993 
2994   Level: advanced
2995 
2996   Note: This is NOT what one uses to update the ghost points before a function evaluation, that should be done at the beginning of your FormFunction()
2997         This is not used by most users.
2998 
2999 .keywords: SNES, update
3000 
3001 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
3002 @*/
3003 PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
3004 {
3005   PetscFunctionBegin;
3006   PetscValidHeaderSpecific(snes, SNES_CLASSID,1);
3007   snes->ops->update = func;
3008   PetscFunctionReturn(0);
3009 }
3010 
3011 #undef __FUNCT__
3012 #define __FUNCT__ "SNESDefaultUpdate"
3013 /*@
3014   SNESDefaultUpdate - The default update function which does nothing.
3015 
3016   Not collective
3017 
3018   Input Parameters:
3019 . snes - The nonlinear solver context
3020 . step - The current step of the iteration
3021 
3022   Level: intermediate
3023 
3024 .keywords: SNES, update
3025 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
3026 @*/
3027 PetscErrorCode  SNESDefaultUpdate(SNES snes, PetscInt step)
3028 {
3029   PetscFunctionBegin;
3030   PetscFunctionReturn(0);
3031 }
3032 
3033 #undef __FUNCT__
3034 #define __FUNCT__ "SNESScaleStep_Private"
3035 /*
3036    SNESScaleStep_Private - Scales a step so that its length is less than the
3037    positive parameter delta.
3038 
3039     Input Parameters:
3040 +   snes - the SNES context
3041 .   y - approximate solution of linear system
3042 .   fnorm - 2-norm of current function
3043 -   delta - trust region size
3044 
3045     Output Parameters:
3046 +   gpnorm - predicted function norm at the new point, assuming local
3047     linearization.  The value is zero if the step lies within the trust
3048     region, and exceeds zero otherwise.
3049 -   ynorm - 2-norm of the step
3050 
3051     Note:
3052     For non-trust region methods such as SNESLS, the parameter delta
3053     is set to be the maximum allowable step size.
3054 
3055 .keywords: SNES, nonlinear, scale, step
3056 */
3057 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
3058 {
3059   PetscReal      nrm;
3060   PetscScalar    cnorm;
3061   PetscErrorCode ierr;
3062 
3063   PetscFunctionBegin;
3064   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3065   PetscValidHeaderSpecific(y,VEC_CLASSID,2);
3066   PetscCheckSameComm(snes,1,y,2);
3067 
3068   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
3069   if (nrm > *delta) {
3070      nrm = *delta/nrm;
3071      *gpnorm = (1.0 - nrm)*(*fnorm);
3072      cnorm = nrm;
3073      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
3074      *ynorm = *delta;
3075   } else {
3076      *gpnorm = 0.0;
3077      *ynorm = nrm;
3078   }
3079   PetscFunctionReturn(0);
3080 }
3081 
3082 #undef __FUNCT__
3083 #define __FUNCT__ "SNESSolve"
3084 /*@C
3085    SNESSolve - Solves a nonlinear system F(x) = b.
3086    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
3087 
3088    Collective on SNES
3089 
3090    Input Parameters:
3091 +  snes - the SNES context
3092 .  b - the constant part of the equation F(x) = b, or PETSC_NULL to use zero.
3093 -  x - the solution vector.
3094 
3095    Notes:
3096    The user should initialize the vector,x, with the initial guess
3097    for the nonlinear solve prior to calling SNESSolve.  In particular,
3098    to employ an initial guess of zero, the user should explicitly set
3099    this vector to zero by calling VecSet().
3100 
3101    Level: beginner
3102 
3103 .keywords: SNES, nonlinear, solve
3104 
3105 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3106 @*/
3107 PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
3108 {
3109   PetscErrorCode ierr;
3110   PetscBool      flg;
3111   char           filename[PETSC_MAX_PATH_LEN];
3112   PetscViewer    viewer;
3113   PetscInt       grid;
3114   Vec            xcreated = PETSC_NULL;
3115 
3116   PetscFunctionBegin;
3117   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3118   if (x) PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3119   if (x) PetscCheckSameComm(snes,1,x,3);
3120   if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2);
3121   if (b) PetscCheckSameComm(snes,1,b,2);
3122 
3123   if (!x && snes->dm) {
3124     ierr = DMCreateGlobalVector(snes->dm,&xcreated);CHKERRQ(ierr);
3125     x    = xcreated;
3126   }
3127 
3128   for (grid=0; grid<snes->gridsequence; grid++) {ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm));CHKERRQ(ierr);}
3129   for (grid=0; grid<snes->gridsequence+1; grid++) {
3130 
3131     /* set solution vector */
3132     if (!grid) {ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);}
3133     ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr);
3134     snes->vec_sol = x;
3135     /* set afine vector if provided */
3136     if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
3137     ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr);
3138     snes->vec_rhs = b;
3139 
3140     ierr = SNESSetUp(snes);CHKERRQ(ierr);
3141 
3142     if (!grid) {
3143       if (snes->ops->computeinitialguess) {
3144         ierr = (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);CHKERRQ(ierr);
3145       } else if (snes->dm) {
3146         PetscBool ig;
3147         ierr = DMHasInitialGuess(snes->dm,&ig);CHKERRQ(ierr);
3148         if (ig) {
3149           ierr = DMComputeInitialGuess(snes->dm,snes->vec_sol);CHKERRQ(ierr);
3150         }
3151       }
3152     }
3153 
3154     if (snes->conv_hist_reset) snes->conv_hist_len = 0;
3155     snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
3156 
3157     ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
3158     ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
3159     ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
3160     if (snes->domainerror){
3161       snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
3162       snes->domainerror = PETSC_FALSE;
3163     }
3164     if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
3165 
3166     ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
3167     if (flg && !PetscPreLoadingOn) {
3168       ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
3169       ierr = SNESView(snes,viewer);CHKERRQ(ierr);
3170       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3171     }
3172 
3173     flg  = PETSC_FALSE;
3174     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr);
3175     if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
3176     if (snes->printreason) {
3177       ierr = PetscViewerASCIIAddTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),((PetscObject)snes)->tablevel);CHKERRQ(ierr);
3178       if (snes->reason > 0) {
3179         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
3180       } else {
3181         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
3182       }
3183       ierr = PetscViewerASCIISubtractTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),((PetscObject)snes)->tablevel);CHKERRQ(ierr);
3184     }
3185 
3186     flg = PETSC_FALSE;
3187     ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view_solution_vtk",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
3188     if (flg) {
3189       PetscViewer viewer;
3190       ierr = PetscViewerCreate(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
3191       ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
3192       ierr = PetscViewerFileSetName(viewer,filename);CHKERRQ(ierr);
3193       ierr = VecView(snes->vec_sol,viewer);CHKERRQ(ierr);
3194       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3195     }
3196 
3197     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3198     if (grid <  snes->gridsequence) {
3199       DM  fine;
3200       Vec xnew;
3201       Mat interp;
3202 
3203       ierr = DMRefine(snes->dm,((PetscObject)snes)->comm,&fine);CHKERRQ(ierr);
3204       ierr = DMCreateInterpolation(snes->dm,fine,&interp,PETSC_NULL);CHKERRQ(ierr);
3205       ierr = DMCreateGlobalVector(fine,&xnew);CHKERRQ(ierr);
3206       ierr = MatInterpolate(interp,x,xnew);CHKERRQ(ierr);
3207       ierr = MatDestroy(&interp);CHKERRQ(ierr);
3208       x    = xnew;
3209 
3210       ierr = SNESReset(snes);CHKERRQ(ierr);
3211       ierr = SNESSetDM(snes,fine);CHKERRQ(ierr);
3212       ierr = DMDestroy(&fine);CHKERRQ(ierr);
3213       ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm));CHKERRQ(ierr);
3214     }
3215   }
3216   ierr = VecDestroy(&xcreated);CHKERRQ(ierr);
3217   PetscFunctionReturn(0);
3218 }
3219 
3220 /* --------- Internal routines for SNES Package --------- */
3221 
3222 #undef __FUNCT__
3223 #define __FUNCT__ "SNESSetType"
3224 /*@C
3225    SNESSetType - Sets the method for the nonlinear solver.
3226 
3227    Collective on SNES
3228 
3229    Input Parameters:
3230 +  snes - the SNES context
3231 -  type - a known method
3232 
3233    Options Database Key:
3234 .  -snes_type <type> - Sets the method; use -help for a list
3235    of available methods (for instance, ls or tr)
3236 
3237    Notes:
3238    See "petsc/include/petscsnes.h" for available methods (for instance)
3239 +    SNESLS - Newton's method with line search
3240      (systems of nonlinear equations)
3241 .    SNESTR - Newton's method with trust region
3242      (systems of nonlinear equations)
3243 
3244   Normally, it is best to use the SNESSetFromOptions() command and then
3245   set the SNES solver type from the options database rather than by using
3246   this routine.  Using the options database provides the user with
3247   maximum flexibility in evaluating the many nonlinear solvers.
3248   The SNESSetType() routine is provided for those situations where it
3249   is necessary to set the nonlinear solver independently of the command
3250   line or options database.  This might be the case, for example, when
3251   the choice of solver changes during the execution of the program,
3252   and the user's application is taking responsibility for choosing the
3253   appropriate method.
3254 
3255   Level: intermediate
3256 
3257 .keywords: SNES, set, type
3258 
3259 .seealso: SNESType, SNESCreate()
3260 
3261 @*/
3262 PetscErrorCode  SNESSetType(SNES snes,const SNESType type)
3263 {
3264   PetscErrorCode ierr,(*r)(SNES);
3265   PetscBool      match;
3266 
3267   PetscFunctionBegin;
3268   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3269   PetscValidCharPointer(type,2);
3270 
3271   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
3272   if (match) PetscFunctionReturn(0);
3273 
3274   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
3275   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
3276   /* Destroy the previous private SNES context */
3277   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
3278   /* Reinitialize function pointers in SNESOps structure */
3279   snes->ops->setup          = 0;
3280   snes->ops->solve          = 0;
3281   snes->ops->view           = 0;
3282   snes->ops->setfromoptions = 0;
3283   snes->ops->destroy        = 0;
3284   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
3285   snes->setupcalled = PETSC_FALSE;
3286   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
3287   ierr = (*r)(snes);CHKERRQ(ierr);
3288 #if defined(PETSC_HAVE_AMS)
3289   if (PetscAMSPublishAll) {
3290     ierr = PetscObjectAMSPublish((PetscObject)snes);CHKERRQ(ierr);
3291   }
3292 #endif
3293   PetscFunctionReturn(0);
3294 }
3295 
3296 
3297 /* --------------------------------------------------------------------- */
3298 #undef __FUNCT__
3299 #define __FUNCT__ "SNESRegisterDestroy"
3300 /*@
3301    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
3302    registered by SNESRegisterDynamic().
3303 
3304    Not Collective
3305 
3306    Level: advanced
3307 
3308 .keywords: SNES, nonlinear, register, destroy
3309 
3310 .seealso: SNESRegisterAll(), SNESRegisterAll()
3311 @*/
3312 PetscErrorCode  SNESRegisterDestroy(void)
3313 {
3314   PetscErrorCode ierr;
3315 
3316   PetscFunctionBegin;
3317   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
3318   SNESRegisterAllCalled = PETSC_FALSE;
3319   PetscFunctionReturn(0);
3320 }
3321 
3322 #undef __FUNCT__
3323 #define __FUNCT__ "SNESGetType"
3324 /*@C
3325    SNESGetType - Gets the SNES method type and name (as a string).
3326 
3327    Not Collective
3328 
3329    Input Parameter:
3330 .  snes - nonlinear solver context
3331 
3332    Output Parameter:
3333 .  type - SNES method (a character string)
3334 
3335    Level: intermediate
3336 
3337 .keywords: SNES, nonlinear, get, type, name
3338 @*/
3339 PetscErrorCode  SNESGetType(SNES snes,const SNESType *type)
3340 {
3341   PetscFunctionBegin;
3342   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3343   PetscValidPointer(type,2);
3344   *type = ((PetscObject)snes)->type_name;
3345   PetscFunctionReturn(0);
3346 }
3347 
3348 #undef __FUNCT__
3349 #define __FUNCT__ "SNESGetSolution"
3350 /*@
3351    SNESGetSolution - Returns the vector where the approximate solution is
3352    stored. This is the fine grid solution when using SNESSetGridSequence().
3353 
3354    Not Collective, but Vec is parallel if SNES is parallel
3355 
3356    Input Parameter:
3357 .  snes - the SNES context
3358 
3359    Output Parameter:
3360 .  x - the solution
3361 
3362    Level: intermediate
3363 
3364 .keywords: SNES, nonlinear, get, solution
3365 
3366 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
3367 @*/
3368 PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
3369 {
3370   PetscFunctionBegin;
3371   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3372   PetscValidPointer(x,2);
3373   *x = snes->vec_sol;
3374   PetscFunctionReturn(0);
3375 }
3376 
3377 #undef __FUNCT__
3378 #define __FUNCT__ "SNESGetSolutionUpdate"
3379 /*@
3380    SNESGetSolutionUpdate - Returns the vector where the solution update is
3381    stored.
3382 
3383    Not Collective, but Vec is parallel if SNES is parallel
3384 
3385    Input Parameter:
3386 .  snes - the SNES context
3387 
3388    Output Parameter:
3389 .  x - the solution update
3390 
3391    Level: advanced
3392 
3393 .keywords: SNES, nonlinear, get, solution, update
3394 
3395 .seealso: SNESGetSolution(), SNESGetFunction()
3396 @*/
3397 PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
3398 {
3399   PetscFunctionBegin;
3400   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3401   PetscValidPointer(x,2);
3402   *x = snes->vec_sol_update;
3403   PetscFunctionReturn(0);
3404 }
3405 
3406 #undef __FUNCT__
3407 #define __FUNCT__ "SNESGetFunction"
3408 /*@C
3409    SNESGetFunction - Returns the vector where the function is stored.
3410 
3411    Not Collective, but Vec is parallel if SNES is parallel. Collective if Vec is requested, but has not been created yet.
3412 
3413    Input Parameter:
3414 .  snes - the SNES context
3415 
3416    Output Parameter:
3417 +  r - the function (or PETSC_NULL)
3418 .  func - the function (or PETSC_NULL)
3419 -  ctx - the function context (or PETSC_NULL)
3420 
3421    Level: advanced
3422 
3423 .keywords: SNES, nonlinear, get, function
3424 
3425 .seealso: SNESSetFunction(), SNESGetSolution()
3426 @*/
3427 PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
3428 {
3429   PetscErrorCode ierr;
3430 
3431   PetscFunctionBegin;
3432   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3433   if (r) {
3434     if (!snes->vec_func) {
3435       if (snes->vec_rhs) {
3436         ierr = VecDuplicate(snes->vec_rhs,&snes->vec_func);CHKERRQ(ierr);
3437       } else if (snes->vec_sol) {
3438         ierr = VecDuplicate(snes->vec_sol,&snes->vec_func);CHKERRQ(ierr);
3439       } else if (snes->dm) {
3440         ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr);
3441       }
3442     }
3443     *r = snes->vec_func;
3444   }
3445   if (func) *func = snes->ops->computefunction;
3446   if (ctx)  *ctx  = snes->funP;
3447   PetscFunctionReturn(0);
3448 }
3449 
3450 /*@C
3451    SNESGetGS - Returns the GS function and context.
3452 
3453    Input Parameter:
3454 .  snes - the SNES context
3455 
3456    Output Parameter:
3457 +  gsfunc - the function (or PETSC_NULL)
3458 -  ctx    - the function context (or PETSC_NULL)
3459 
3460    Level: advanced
3461 
3462 .keywords: SNES, nonlinear, get, function
3463 
3464 .seealso: SNESSetGS(), SNESGetFunction()
3465 @*/
3466 
3467 #undef __FUNCT__
3468 #define __FUNCT__ "SNESGetGS"
3469 PetscErrorCode SNESGetGS (SNES snes, PetscErrorCode(**func)(SNES, Vec, Vec, void*), void ** ctx)
3470 {
3471   PetscFunctionBegin;
3472   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3473   if (func) *func = snes->ops->computegs;
3474   if (ctx)  *ctx  = snes->funP;
3475   PetscFunctionReturn(0);
3476 }
3477 
3478 #undef __FUNCT__
3479 #define __FUNCT__ "SNESSetOptionsPrefix"
3480 /*@C
3481    SNESSetOptionsPrefix - Sets the prefix used for searching for all
3482    SNES options in the database.
3483 
3484    Logically Collective on SNES
3485 
3486    Input Parameter:
3487 +  snes - the SNES context
3488 -  prefix - the prefix to prepend to all option names
3489 
3490    Notes:
3491    A hyphen (-) must NOT be given at the beginning of the prefix name.
3492    The first character of all runtime options is AUTOMATICALLY the hyphen.
3493 
3494    Level: advanced
3495 
3496 .keywords: SNES, set, options, prefix, database
3497 
3498 .seealso: SNESSetFromOptions()
3499 @*/
3500 PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
3501 {
3502   PetscErrorCode ierr;
3503 
3504   PetscFunctionBegin;
3505   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3506   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
3507   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
3508   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
3509   PetscFunctionReturn(0);
3510 }
3511 
3512 #undef __FUNCT__
3513 #define __FUNCT__ "SNESAppendOptionsPrefix"
3514 /*@C
3515    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
3516    SNES options in the database.
3517 
3518    Logically Collective on SNES
3519 
3520    Input Parameters:
3521 +  snes - the SNES context
3522 -  prefix - the prefix to prepend to all option names
3523 
3524    Notes:
3525    A hyphen (-) must NOT be given at the beginning of the prefix name.
3526    The first character of all runtime options is AUTOMATICALLY the hyphen.
3527 
3528    Level: advanced
3529 
3530 .keywords: SNES, append, options, prefix, database
3531 
3532 .seealso: SNESGetOptionsPrefix()
3533 @*/
3534 PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
3535 {
3536   PetscErrorCode ierr;
3537 
3538   PetscFunctionBegin;
3539   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3540   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
3541   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
3542   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
3543   PetscFunctionReturn(0);
3544 }
3545 
3546 #undef __FUNCT__
3547 #define __FUNCT__ "SNESGetOptionsPrefix"
3548 /*@C
3549    SNESGetOptionsPrefix - Sets the prefix used for searching for all
3550    SNES options in the database.
3551 
3552    Not Collective
3553 
3554    Input Parameter:
3555 .  snes - the SNES context
3556 
3557    Output Parameter:
3558 .  prefix - pointer to the prefix string used
3559 
3560    Notes: On the fortran side, the user should pass in a string 'prefix' of
3561    sufficient length to hold the prefix.
3562 
3563    Level: advanced
3564 
3565 .keywords: SNES, get, options, prefix, database
3566 
3567 .seealso: SNESAppendOptionsPrefix()
3568 @*/
3569 PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
3570 {
3571   PetscErrorCode ierr;
3572 
3573   PetscFunctionBegin;
3574   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3575   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
3576   PetscFunctionReturn(0);
3577 }
3578 
3579 
3580 #undef __FUNCT__
3581 #define __FUNCT__ "SNESRegister"
3582 /*@C
3583   SNESRegister - See SNESRegisterDynamic()
3584 
3585   Level: advanced
3586 @*/
3587 PetscErrorCode  SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
3588 {
3589   char           fullname[PETSC_MAX_PATH_LEN];
3590   PetscErrorCode ierr;
3591 
3592   PetscFunctionBegin;
3593   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
3594   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
3595   PetscFunctionReturn(0);
3596 }
3597 
3598 #undef __FUNCT__
3599 #define __FUNCT__ "SNESTestLocalMin"
3600 PetscErrorCode  SNESTestLocalMin(SNES snes)
3601 {
3602   PetscErrorCode ierr;
3603   PetscInt       N,i,j;
3604   Vec            u,uh,fh;
3605   PetscScalar    value;
3606   PetscReal      norm;
3607 
3608   PetscFunctionBegin;
3609   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
3610   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
3611   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
3612 
3613   /* currently only works for sequential */
3614   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
3615   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
3616   for (i=0; i<N; i++) {
3617     ierr = VecCopy(u,uh);CHKERRQ(ierr);
3618     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
3619     for (j=-10; j<11; j++) {
3620       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
3621       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
3622       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
3623       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
3624       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
3625       value = -value;
3626       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
3627     }
3628   }
3629   ierr = VecDestroy(&uh);CHKERRQ(ierr);
3630   ierr = VecDestroy(&fh);CHKERRQ(ierr);
3631   PetscFunctionReturn(0);
3632 }
3633 
3634 #undef __FUNCT__
3635 #define __FUNCT__ "SNESKSPSetUseEW"
3636 /*@
3637    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
3638    computing relative tolerance for linear solvers within an inexact
3639    Newton method.
3640 
3641    Logically Collective on SNES
3642 
3643    Input Parameters:
3644 +  snes - SNES context
3645 -  flag - PETSC_TRUE or PETSC_FALSE
3646 
3647     Options Database:
3648 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
3649 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
3650 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
3651 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
3652 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
3653 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
3654 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
3655 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
3656 
3657    Notes:
3658    Currently, the default is to use a constant relative tolerance for
3659    the inner linear solvers.  Alternatively, one can use the
3660    Eisenstat-Walker method, where the relative convergence tolerance
3661    is reset at each Newton iteration according progress of the nonlinear
3662    solver.
3663 
3664    Level: advanced
3665 
3666    Reference:
3667    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
3668    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
3669 
3670 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
3671 
3672 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
3673 @*/
3674 PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool  flag)
3675 {
3676   PetscFunctionBegin;
3677   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3678   PetscValidLogicalCollectiveBool(snes,flag,2);
3679   snes->ksp_ewconv = flag;
3680   PetscFunctionReturn(0);
3681 }
3682 
3683 #undef __FUNCT__
3684 #define __FUNCT__ "SNESKSPGetUseEW"
3685 /*@
3686    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
3687    for computing relative tolerance for linear solvers within an
3688    inexact Newton method.
3689 
3690    Not Collective
3691 
3692    Input Parameter:
3693 .  snes - SNES context
3694 
3695    Output Parameter:
3696 .  flag - PETSC_TRUE or PETSC_FALSE
3697 
3698    Notes:
3699    Currently, the default is to use a constant relative tolerance for
3700    the inner linear solvers.  Alternatively, one can use the
3701    Eisenstat-Walker method, where the relative convergence tolerance
3702    is reset at each Newton iteration according progress of the nonlinear
3703    solver.
3704 
3705    Level: advanced
3706 
3707    Reference:
3708    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
3709    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
3710 
3711 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
3712 
3713 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
3714 @*/
3715 PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
3716 {
3717   PetscFunctionBegin;
3718   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3719   PetscValidPointer(flag,2);
3720   *flag = snes->ksp_ewconv;
3721   PetscFunctionReturn(0);
3722 }
3723 
3724 #undef __FUNCT__
3725 #define __FUNCT__ "SNESKSPSetParametersEW"
3726 /*@
3727    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
3728    convergence criteria for the linear solvers within an inexact
3729    Newton method.
3730 
3731    Logically Collective on SNES
3732 
3733    Input Parameters:
3734 +    snes - SNES context
3735 .    version - version 1, 2 (default is 2) or 3
3736 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
3737 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
3738 .    gamma - multiplicative factor for version 2 rtol computation
3739              (0 <= gamma2 <= 1)
3740 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
3741 .    alpha2 - power for safeguard
3742 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
3743 
3744    Note:
3745    Version 3 was contributed by Luis Chacon, June 2006.
3746 
3747    Use PETSC_DEFAULT to retain the default for any of the parameters.
3748 
3749    Level: advanced
3750 
3751    Reference:
3752    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
3753    inexact Newton method", Utah State University Math. Stat. Dept. Res.
3754    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
3755 
3756 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
3757 
3758 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
3759 @*/
3760 PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
3761 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
3762 {
3763   SNESKSPEW *kctx;
3764   PetscFunctionBegin;
3765   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3766   kctx = (SNESKSPEW*)snes->kspconvctx;
3767   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
3768   PetscValidLogicalCollectiveInt(snes,version,2);
3769   PetscValidLogicalCollectiveReal(snes,rtol_0,3);
3770   PetscValidLogicalCollectiveReal(snes,rtol_max,4);
3771   PetscValidLogicalCollectiveReal(snes,gamma,5);
3772   PetscValidLogicalCollectiveReal(snes,alpha,6);
3773   PetscValidLogicalCollectiveReal(snes,alpha2,7);
3774   PetscValidLogicalCollectiveReal(snes,threshold,8);
3775 
3776   if (version != PETSC_DEFAULT)   kctx->version   = version;
3777   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
3778   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
3779   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
3780   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
3781   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
3782   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
3783 
3784   if (kctx->version < 1 || kctx->version > 3) {
3785     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
3786   }
3787   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
3788     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
3789   }
3790   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
3791     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
3792   }
3793   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
3794     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
3795   }
3796   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
3797     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
3798   }
3799   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
3800     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
3801   }
3802   PetscFunctionReturn(0);
3803 }
3804 
3805 #undef __FUNCT__
3806 #define __FUNCT__ "SNESKSPGetParametersEW"
3807 /*@
3808    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
3809    convergence criteria for the linear solvers within an inexact
3810    Newton method.
3811 
3812    Not Collective
3813 
3814    Input Parameters:
3815      snes - SNES context
3816 
3817    Output Parameters:
3818 +    version - version 1, 2 (default is 2) or 3
3819 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
3820 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
3821 .    gamma - multiplicative factor for version 2 rtol computation
3822              (0 <= gamma2 <= 1)
3823 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
3824 .    alpha2 - power for safeguard
3825 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
3826 
3827    Level: advanced
3828 
3829 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
3830 
3831 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
3832 @*/
3833 PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
3834 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
3835 {
3836   SNESKSPEW *kctx;
3837   PetscFunctionBegin;
3838   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3839   kctx = (SNESKSPEW*)snes->kspconvctx;
3840   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
3841   if(version)   *version   = kctx->version;
3842   if(rtol_0)    *rtol_0    = kctx->rtol_0;
3843   if(rtol_max)  *rtol_max  = kctx->rtol_max;
3844   if(gamma)     *gamma     = kctx->gamma;
3845   if(alpha)     *alpha     = kctx->alpha;
3846   if(alpha2)    *alpha2    = kctx->alpha2;
3847   if(threshold) *threshold = kctx->threshold;
3848   PetscFunctionReturn(0);
3849 }
3850 
3851 #undef __FUNCT__
3852 #define __FUNCT__ "SNESKSPEW_PreSolve"
3853 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
3854 {
3855   PetscErrorCode ierr;
3856   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
3857   PetscReal      rtol=PETSC_DEFAULT,stol;
3858 
3859   PetscFunctionBegin;
3860   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
3861   if (!snes->iter) { /* first time in, so use the original user rtol */
3862     rtol = kctx->rtol_0;
3863   } else {
3864     if (kctx->version == 1) {
3865       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
3866       if (rtol < 0.0) rtol = -rtol;
3867       stol = pow(kctx->rtol_last,kctx->alpha2);
3868       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
3869     } else if (kctx->version == 2) {
3870       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
3871       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
3872       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
3873     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
3874       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
3875       /* safeguard: avoid sharp decrease of rtol */
3876       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
3877       stol = PetscMax(rtol,stol);
3878       rtol = PetscMin(kctx->rtol_0,stol);
3879       /* safeguard: avoid oversolving */
3880       stol = kctx->gamma*(snes->ttol)/snes->norm;
3881       stol = PetscMax(rtol,stol);
3882       rtol = PetscMin(kctx->rtol_0,stol);
3883     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
3884   }
3885   /* safeguard: avoid rtol greater than one */
3886   rtol = PetscMin(rtol,kctx->rtol_max);
3887   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
3888   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
3889   PetscFunctionReturn(0);
3890 }
3891 
3892 #undef __FUNCT__
3893 #define __FUNCT__ "SNESKSPEW_PostSolve"
3894 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
3895 {
3896   PetscErrorCode ierr;
3897   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
3898   PCSide         pcside;
3899   Vec            lres;
3900 
3901   PetscFunctionBegin;
3902   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
3903   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
3904   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
3905   if (kctx->version == 1) {
3906     ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr);
3907     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
3908       /* KSP residual is true linear residual */
3909       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
3910     } else {
3911       /* KSP residual is preconditioned residual */
3912       /* compute true linear residual norm */
3913       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
3914       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
3915       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
3916       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
3917       ierr = VecDestroy(&lres);CHKERRQ(ierr);
3918     }
3919   }
3920   PetscFunctionReturn(0);
3921 }
3922 
3923 #undef __FUNCT__
3924 #define __FUNCT__ "SNES_KSPSolve"
3925 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
3926 {
3927   PetscErrorCode ierr;
3928 
3929   PetscFunctionBegin;
3930   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
3931   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
3932   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
3933   PetscFunctionReturn(0);
3934 }
3935 
3936 #undef __FUNCT__
3937 #define __FUNCT__ "SNESSetDM"
3938 /*@
3939    SNESSetDM - Sets the DM that may be used by some preconditioners
3940 
3941    Logically Collective on SNES
3942 
3943    Input Parameters:
3944 +  snes - the preconditioner context
3945 -  dm - the dm
3946 
3947    Level: intermediate
3948 
3949 
3950 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
3951 @*/
3952 PetscErrorCode  SNESSetDM(SNES snes,DM dm)
3953 {
3954   PetscErrorCode ierr;
3955   KSP            ksp;
3956 
3957   PetscFunctionBegin;
3958   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3959   if (dm) {ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);}
3960   ierr = DMDestroy(&snes->dm);CHKERRQ(ierr);
3961   snes->dm = dm;
3962   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
3963   ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr);
3964   ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr);
3965   if (snes->pc) {
3966     ierr = SNESSetDM(snes->pc, snes->dm);CHKERRQ(ierr);
3967   }
3968   PetscFunctionReturn(0);
3969 }
3970 
3971 #undef __FUNCT__
3972 #define __FUNCT__ "SNESGetDM"
3973 /*@
3974    SNESGetDM - Gets the DM that may be used by some preconditioners
3975 
3976    Not Collective but DM obtained is parallel on SNES
3977 
3978    Input Parameter:
3979 . snes - the preconditioner context
3980 
3981    Output Parameter:
3982 .  dm - the dm
3983 
3984    Level: intermediate
3985 
3986 
3987 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
3988 @*/
3989 PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
3990 {
3991   PetscFunctionBegin;
3992   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3993   *dm = snes->dm;
3994   PetscFunctionReturn(0);
3995 }
3996 
3997 #undef __FUNCT__
3998 #define __FUNCT__ "SNESSetPC"
3999 /*@
4000   SNESSetPC - Sets the nonlinear preconditioner to be used.
4001 
4002   Collective on SNES
4003 
4004   Input Parameters:
4005 + snes - iterative context obtained from SNESCreate()
4006 - pc   - the preconditioner object
4007 
4008   Notes:
4009   Use SNESGetPC() to retrieve the preconditioner context (for example,
4010   to configure it using the API).
4011 
4012   Level: developer
4013 
4014 .keywords: SNES, set, precondition
4015 .seealso: SNESGetPC()
4016 @*/
4017 PetscErrorCode SNESSetPC(SNES snes, SNES pc)
4018 {
4019   PetscErrorCode ierr;
4020 
4021   PetscFunctionBegin;
4022   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4023   PetscValidHeaderSpecific(pc, SNES_CLASSID, 2);
4024   PetscCheckSameComm(snes, 1, pc, 2);
4025   ierr = PetscObjectReference((PetscObject) pc);CHKERRQ(ierr);
4026   ierr = SNESDestroy(&snes->pc);CHKERRQ(ierr);
4027   snes->pc = pc;
4028   ierr = PetscLogObjectParent(snes, snes->pc);CHKERRQ(ierr);
4029   PetscFunctionReturn(0);
4030 }
4031 
4032 #undef __FUNCT__
4033 #define __FUNCT__ "SNESGetPC"
4034 /*@
4035   SNESGetPC - Returns a pointer to the nonlinear preconditioning context set with SNESSetPC().
4036 
4037   Not Collective
4038 
4039   Input Parameter:
4040 . snes - iterative context obtained from SNESCreate()
4041 
4042   Output Parameter:
4043 . pc - preconditioner context
4044 
4045   Level: developer
4046 
4047 .keywords: SNES, get, preconditioner
4048 .seealso: SNESSetPC()
4049 @*/
4050 PetscErrorCode SNESGetPC(SNES snes, SNES *pc)
4051 {
4052   PetscErrorCode ierr;
4053 
4054   PetscFunctionBegin;
4055   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
4056   PetscValidPointer(pc, 2);
4057   if (!snes->pc) {
4058     ierr = SNESCreate(((PetscObject) snes)->comm, &snes->pc);CHKERRQ(ierr);
4059     ierr = PetscObjectIncrementTabLevel((PetscObject) snes->pc, (PetscObject) snes, 1);CHKERRQ(ierr);
4060     ierr = PetscLogObjectParent(snes, snes->pc);CHKERRQ(ierr);
4061   }
4062   *pc = snes->pc;
4063   PetscFunctionReturn(0);
4064 }
4065 
4066 #if defined(PETSC_HAVE_MATLAB_ENGINE)
4067 #include <mex.h>
4068 
4069 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
4070 
4071 #undef __FUNCT__
4072 #define __FUNCT__ "SNESComputeFunction_Matlab"
4073 /*
4074    SNESComputeFunction_Matlab - Calls the function that has been set with
4075                          SNESSetFunctionMatlab().
4076 
4077    Collective on SNES
4078 
4079    Input Parameters:
4080 +  snes - the SNES context
4081 -  x - input vector
4082 
4083    Output Parameter:
4084 .  y - function vector, as set by SNESSetFunction()
4085 
4086    Notes:
4087    SNESComputeFunction() is typically used within nonlinear solvers
4088    implementations, so most users would not generally call this routine
4089    themselves.
4090 
4091    Level: developer
4092 
4093 .keywords: SNES, nonlinear, compute, function
4094 
4095 .seealso: SNESSetFunction(), SNESGetFunction()
4096 */
4097 PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
4098 {
4099   PetscErrorCode    ierr;
4100   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
4101   int               nlhs = 1,nrhs = 5;
4102   mxArray	    *plhs[1],*prhs[5];
4103   long long int     lx = 0,ly = 0,ls = 0;
4104 
4105   PetscFunctionBegin;
4106   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4107   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
4108   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
4109   PetscCheckSameComm(snes,1,x,2);
4110   PetscCheckSameComm(snes,1,y,3);
4111 
4112   /* call Matlab function in ctx with arguments x and y */
4113 
4114   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
4115   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
4116   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
4117   prhs[0] =  mxCreateDoubleScalar((double)ls);
4118   prhs[1] =  mxCreateDoubleScalar((double)lx);
4119   prhs[2] =  mxCreateDoubleScalar((double)ly);
4120   prhs[3] =  mxCreateString(sctx->funcname);
4121   prhs[4] =  sctx->ctx;
4122   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");CHKERRQ(ierr);
4123   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4124   mxDestroyArray(prhs[0]);
4125   mxDestroyArray(prhs[1]);
4126   mxDestroyArray(prhs[2]);
4127   mxDestroyArray(prhs[3]);
4128   mxDestroyArray(plhs[0]);
4129   PetscFunctionReturn(0);
4130 }
4131 
4132 
4133 #undef __FUNCT__
4134 #define __FUNCT__ "SNESSetFunctionMatlab"
4135 /*
4136    SNESSetFunctionMatlab - Sets the function evaluation routine and function
4137    vector for use by the SNES routines in solving systems of nonlinear
4138    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
4139 
4140    Logically Collective on SNES
4141 
4142    Input Parameters:
4143 +  snes - the SNES context
4144 .  r - vector to store function value
4145 -  func - function evaluation routine
4146 
4147    Calling sequence of func:
4148 $    func (SNES snes,Vec x,Vec f,void *ctx);
4149 
4150 
4151    Notes:
4152    The Newton-like methods typically solve linear systems of the form
4153 $      f'(x) x = -f(x),
4154    where f'(x) denotes the Jacobian matrix and f(x) is the function.
4155 
4156    Level: beginner
4157 
4158 .keywords: SNES, nonlinear, set, function
4159 
4160 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4161 */
4162 PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *func,mxArray *ctx)
4163 {
4164   PetscErrorCode    ierr;
4165   SNESMatlabContext *sctx;
4166 
4167   PetscFunctionBegin;
4168   /* currently sctx is memory bleed */
4169   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
4170   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4171   /*
4172      This should work, but it doesn't
4173   sctx->ctx = ctx;
4174   mexMakeArrayPersistent(sctx->ctx);
4175   */
4176   sctx->ctx = mxDuplicateArray(ctx);
4177   ierr = SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);CHKERRQ(ierr);
4178   PetscFunctionReturn(0);
4179 }
4180 
4181 #undef __FUNCT__
4182 #define __FUNCT__ "SNESComputeJacobian_Matlab"
4183 /*
4184    SNESComputeJacobian_Matlab - Calls the function that has been set with
4185                          SNESSetJacobianMatlab().
4186 
4187    Collective on SNES
4188 
4189    Input Parameters:
4190 +  snes - the SNES context
4191 .  x - input vector
4192 .  A, B - the matrices
4193 -  ctx - user context
4194 
4195    Output Parameter:
4196 .  flag - structure of the matrix
4197 
4198    Level: developer
4199 
4200 .keywords: SNES, nonlinear, compute, function
4201 
4202 .seealso: SNESSetFunction(), SNESGetFunction()
4203 @*/
4204 PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag, void *ctx)
4205 {
4206   PetscErrorCode    ierr;
4207   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
4208   int               nlhs = 2,nrhs = 6;
4209   mxArray	    *plhs[2],*prhs[6];
4210   long long int     lx = 0,lA = 0,ls = 0, lB = 0;
4211 
4212   PetscFunctionBegin;
4213   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4214   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
4215 
4216   /* call Matlab function in ctx with arguments x and y */
4217 
4218   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
4219   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
4220   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
4221   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
4222   prhs[0] =  mxCreateDoubleScalar((double)ls);
4223   prhs[1] =  mxCreateDoubleScalar((double)lx);
4224   prhs[2] =  mxCreateDoubleScalar((double)lA);
4225   prhs[3] =  mxCreateDoubleScalar((double)lB);
4226   prhs[4] =  mxCreateString(sctx->funcname);
4227   prhs[5] =  sctx->ctx;
4228   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");CHKERRQ(ierr);
4229   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4230   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
4231   mxDestroyArray(prhs[0]);
4232   mxDestroyArray(prhs[1]);
4233   mxDestroyArray(prhs[2]);
4234   mxDestroyArray(prhs[3]);
4235   mxDestroyArray(prhs[4]);
4236   mxDestroyArray(plhs[0]);
4237   mxDestroyArray(plhs[1]);
4238   PetscFunctionReturn(0);
4239 }
4240 
4241 
4242 #undef __FUNCT__
4243 #define __FUNCT__ "SNESSetJacobianMatlab"
4244 /*
4245    SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
4246    vector for use by the SNES routines in solving systems of nonlinear
4247    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
4248 
4249    Logically Collective on SNES
4250 
4251    Input Parameters:
4252 +  snes - the SNES context
4253 .  A,B - Jacobian matrices
4254 .  func - function evaluation routine
4255 -  ctx - user context
4256 
4257    Calling sequence of func:
4258 $    flag = func (SNES snes,Vec x,Mat A,Mat B,void *ctx);
4259 
4260 
4261    Level: developer
4262 
4263 .keywords: SNES, nonlinear, set, function
4264 
4265 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4266 */
4267 PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *func,mxArray *ctx)
4268 {
4269   PetscErrorCode    ierr;
4270   SNESMatlabContext *sctx;
4271 
4272   PetscFunctionBegin;
4273   /* currently sctx is memory bleed */
4274   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
4275   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4276   /*
4277      This should work, but it doesn't
4278   sctx->ctx = ctx;
4279   mexMakeArrayPersistent(sctx->ctx);
4280   */
4281   sctx->ctx = mxDuplicateArray(ctx);
4282   ierr = SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
4283   PetscFunctionReturn(0);
4284 }
4285 
4286 #undef __FUNCT__
4287 #define __FUNCT__ "SNESMonitor_Matlab"
4288 /*
4289    SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
4290 
4291    Collective on SNES
4292 
4293 .seealso: SNESSetFunction(), SNESGetFunction()
4294 @*/
4295 PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
4296 {
4297   PetscErrorCode  ierr;
4298   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
4299   int             nlhs = 1,nrhs = 6;
4300   mxArray	  *plhs[1],*prhs[6];
4301   long long int   lx = 0,ls = 0;
4302   Vec             x=snes->vec_sol;
4303 
4304   PetscFunctionBegin;
4305   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4306 
4307   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
4308   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
4309   prhs[0] =  mxCreateDoubleScalar((double)ls);
4310   prhs[1] =  mxCreateDoubleScalar((double)it);
4311   prhs[2] =  mxCreateDoubleScalar((double)fnorm);
4312   prhs[3] =  mxCreateDoubleScalar((double)lx);
4313   prhs[4] =  mxCreateString(sctx->funcname);
4314   prhs[5] =  sctx->ctx;
4315   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");CHKERRQ(ierr);
4316   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4317   mxDestroyArray(prhs[0]);
4318   mxDestroyArray(prhs[1]);
4319   mxDestroyArray(prhs[2]);
4320   mxDestroyArray(prhs[3]);
4321   mxDestroyArray(prhs[4]);
4322   mxDestroyArray(plhs[0]);
4323   PetscFunctionReturn(0);
4324 }
4325 
4326 
4327 #undef __FUNCT__
4328 #define __FUNCT__ "SNESMonitorSetMatlab"
4329 /*
4330    SNESMonitorSetMatlab - Sets the monitor function from MATLAB
4331 
4332    Level: developer
4333 
4334 .keywords: SNES, nonlinear, set, function
4335 
4336 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4337 */
4338 PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *func,mxArray *ctx)
4339 {
4340   PetscErrorCode    ierr;
4341   SNESMatlabContext *sctx;
4342 
4343   PetscFunctionBegin;
4344   /* currently sctx is memory bleed */
4345   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
4346   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4347   /*
4348      This should work, but it doesn't
4349   sctx->ctx = ctx;
4350   mexMakeArrayPersistent(sctx->ctx);
4351   */
4352   sctx->ctx = mxDuplicateArray(ctx);
4353   ierr = SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
4354   PetscFunctionReturn(0);
4355 }
4356 
4357 #endif
4358