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