xref: /petsc/src/snes/interface/snes.c (revision 701cf23dd8404bcbcedcc00b2ad471028c8837a0)
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
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.
1374 
1375    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
1376    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
1377    different please contact us at petsc-dev@mcs.anl.gov and we'll have an entirely new argument :-).
1378 
1379    Level: beginner
1380 
1381 .keywords: SNES, nonlinear, set, function
1382 
1383 .seealso: SNESGetFunction(), SNESSetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetPicard()
1384 @*/
1385 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)
1386 {
1387   PetscErrorCode ierr;
1388   PetscFunctionBegin;
1389   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1390   snes->ops->computepfunction = func;
1391   snes->ops->computepjacobian = mfunc;
1392   ierr = SNESSetFunction(snes,r,SNESPicardComputeFunction,ctx);CHKERRQ(ierr);
1393   ierr = SNESSetJacobian(snes,jmat,mat,SNESPicardComputeJacobian,ctx);CHKERRQ(ierr);
1394   PetscFunctionReturn(0);
1395 }
1396 
1397 #undef __FUNCT__
1398 #define __FUNCT__ "SNESSetComputeInitialGuess"
1399 /*@C
1400    SNESSetComputeInitialGuess - Sets a routine used to compute an initial guess for the problem
1401 
1402    Logically Collective on SNES
1403 
1404    Input Parameters:
1405 +  snes - the SNES context
1406 .  func - function evaluation routine
1407 -  ctx - [optional] user-defined context for private data for the
1408          function evaluation routine (may be PETSC_NULL)
1409 
1410    Calling sequence of func:
1411 $    func (SNES snes,Vec x,void *ctx);
1412 
1413 .  f - function vector
1414 -  ctx - optional user-defined function context
1415 
1416    Level: intermediate
1417 
1418 .keywords: SNES, nonlinear, set, function
1419 
1420 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
1421 @*/
1422 PetscErrorCode  SNESSetComputeInitialGuess(SNES snes,PetscErrorCode (*func)(SNES,Vec,void*),void *ctx)
1423 {
1424   PetscFunctionBegin;
1425   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1426   if (func) snes->ops->computeinitialguess = func;
1427   if (ctx)  snes->initialguessP            = ctx;
1428   PetscFunctionReturn(0);
1429 }
1430 
1431 /* --------------------------------------------------------------- */
1432 #undef __FUNCT__
1433 #define __FUNCT__ "SNESGetRhs"
1434 /*@C
1435    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
1436    it assumes a zero right hand side.
1437 
1438    Logically Collective on SNES
1439 
1440    Input Parameter:
1441 .  snes - the SNES context
1442 
1443    Output Parameter:
1444 .  rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null
1445 
1446    Level: intermediate
1447 
1448 .keywords: SNES, nonlinear, get, function, right hand side
1449 
1450 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
1451 @*/
1452 PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
1453 {
1454   PetscFunctionBegin;
1455   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1456   PetscValidPointer(rhs,2);
1457   *rhs = snes->vec_rhs;
1458   PetscFunctionReturn(0);
1459 }
1460 
1461 #undef __FUNCT__
1462 #define __FUNCT__ "SNESComputeFunction"
1463 /*@
1464    SNESComputeFunction - Calls the function that has been set with
1465                          SNESSetFunction().
1466 
1467    Collective on SNES
1468 
1469    Input Parameters:
1470 +  snes - the SNES context
1471 -  x - input vector
1472 
1473    Output Parameter:
1474 .  y - function vector, as set by SNESSetFunction()
1475 
1476    Notes:
1477    SNESComputeFunction() is typically used within nonlinear solvers
1478    implementations, so most users would not generally call this routine
1479    themselves.
1480 
1481    Level: developer
1482 
1483 .keywords: SNES, nonlinear, compute, function
1484 
1485 .seealso: SNESSetFunction(), SNESGetFunction()
1486 @*/
1487 PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
1488 {
1489   PetscErrorCode ierr;
1490 
1491   PetscFunctionBegin;
1492   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1493   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1494   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1495   PetscCheckSameComm(snes,1,x,2);
1496   PetscCheckSameComm(snes,1,y,3);
1497 
1498   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1499   if (snes->ops->computefunction) {
1500     PetscStackPush("SNES user function");
1501     ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
1502     PetscStackPop;
1503   } else if (snes->dm) {
1504     ierr = DMComputeFunction(snes->dm,x,y);CHKERRQ(ierr);
1505   } else if (snes->vec_rhs) {
1506     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
1507   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() or SNESSetDM() before SNESComputeFunction(), likely called from SNESSolve().");
1508   if (snes->vec_rhs) {
1509     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
1510   }
1511   snes->nfuncs++;
1512   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1513   PetscFunctionReturn(0);
1514 }
1515 
1516 #undef __FUNCT__
1517 #define __FUNCT__ "SNESComputeGS"
1518 /*@
1519    SNESComputeGS - Calls the Gauss-Seidel function that has been set with
1520                    SNESSetGS().
1521 
1522    Collective on SNES
1523 
1524    Input Parameters:
1525 +  snes - the SNES context
1526 .  x - input vector
1527 -  b - rhs vector
1528 
1529    Output Parameter:
1530 .  x - new solution vector
1531 
1532    Notes:
1533    SNESComputeGS() is typically used within composed nonlinear solver
1534    implementations, so most users would not generally call this routine
1535    themselves.
1536 
1537    Level: developer
1538 
1539 .keywords: SNES, nonlinear, compute, function
1540 
1541 .seealso: SNESSetGS(), SNESComputeFunction()
1542 @*/
1543 PetscErrorCode  SNESComputeGS(SNES snes,Vec b,Vec x)
1544 {
1545   PetscErrorCode ierr;
1546 
1547   PetscFunctionBegin;
1548   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1549   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
1550   if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,3);
1551   PetscCheckSameComm(snes,1,x,2);
1552   if(b) PetscCheckSameComm(snes,1,b,3);
1553   ierr = PetscLogEventBegin(SNES_GSEval,snes,x,b,0);CHKERRQ(ierr);
1554   if (snes->ops->computegs) {
1555     PetscStackPush("SNES user GS");
1556     ierr = (*snes->ops->computegs)(snes,x,b,snes->gsP);CHKERRQ(ierr);
1557     PetscStackPop;
1558   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetGS() before SNESComputeGS(), likely called from SNESSolve().");
1559   ierr = PetscLogEventEnd(SNES_GSEval,snes,x,b,0);CHKERRQ(ierr);
1560   PetscFunctionReturn(0);
1561 }
1562 
1563 
1564 #undef __FUNCT__
1565 #define __FUNCT__ "SNESComputeJacobian"
1566 /*@
1567    SNESComputeJacobian - Computes the Jacobian matrix that has been
1568    set with SNESSetJacobian().
1569 
1570    Collective on SNES and Mat
1571 
1572    Input Parameters:
1573 +  snes - the SNES context
1574 -  x - input vector
1575 
1576    Output Parameters:
1577 +  A - Jacobian matrix
1578 .  B - optional preconditioning matrix
1579 -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1580 
1581   Options Database Keys:
1582 +    -snes_lag_preconditioner <lag>
1583 .    -snes_lag_jacobian <lag>
1584 .    -snes_compare_explicit - Compare the computed Jacobian to the finite difference Jacobian and output the differences
1585 .    -snes_compare_explicit_draw  - Compare the computed Jacobian to the finite difference Jacobian and draw the result
1586 .    -snes_compare_explicit_contour  - Compare the computed Jacobian to the finite difference Jacobian and draw a contour plot with the result
1587 .    -snes_compare_operator  - Make the comparison options above use the operator instead of the preconditioning matrix
1588 .    -snes_compare_coloring - Compute the finite differece Jacobian using coloring and display norms of difference
1589 .    -snes_compare_coloring_display - Compute the finite differece Jacobian using coloring and display verbose differences
1590 .    -snes_compare_coloring_threshold - Display only those matrix entries that differ by more than a given threshold
1591 .    -snes_compare_coloring_threshold_atol - Absolute tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
1592 .    -snes_compare_coloring_threshold_rtol - Relative tolerance for difference in matrix entries to be displayed by -snes_compare_coloring_threshold
1593 .    -snes_compare_coloring_draw - Compute the finite differece Jacobian using coloring and draw differences
1594 -    -snes_compare_coloring_draw_contour - Compute the finite differece Jacobian using coloring and show contours of matrices and differences
1595 
1596 
1597    Notes:
1598    Most users should not need to explicitly call this routine, as it
1599    is used internally within the nonlinear solvers.
1600 
1601    See KSPSetOperators() for important information about setting the
1602    flag parameter.
1603 
1604    Level: developer
1605 
1606 .keywords: SNES, compute, Jacobian, matrix
1607 
1608 .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
1609 @*/
1610 PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
1611 {
1612   PetscErrorCode ierr;
1613   PetscBool      flag;
1614 
1615   PetscFunctionBegin;
1616   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1617   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1618   PetscValidPointer(flg,5);
1619   PetscCheckSameComm(snes,1,X,2);
1620   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1621 
1622   /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
1623 
1624   if (snes->lagjacobian == -2) {
1625     snes->lagjacobian = -1;
1626     ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr);
1627   } else if (snes->lagjacobian == -1) {
1628     *flg = SAME_PRECONDITIONER;
1629     ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr);
1630     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1631     if (flag) {
1632       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1633       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1634     }
1635     PetscFunctionReturn(0);
1636   } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) {
1637     *flg = SAME_PRECONDITIONER;
1638     ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr);
1639     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1640     if (flag) {
1641       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1642       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1643     }
1644     PetscFunctionReturn(0);
1645   }
1646 
1647   *flg = DIFFERENT_NONZERO_PATTERN;
1648   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1649   PetscStackPush("SNES user Jacobian function");
1650   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1651   PetscStackPop;
1652   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1653 
1654   if (snes->lagpreconditioner == -2) {
1655     ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr);
1656     snes->lagpreconditioner = -1;
1657   } else if (snes->lagpreconditioner == -1) {
1658     *flg = SAME_PRECONDITIONER;
1659     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1660   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1661     *flg = SAME_PRECONDITIONER;
1662     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1663   }
1664 
1665   /* make sure user returned a correct Jacobian and preconditioner */
1666   /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
1667     PetscValidHeaderSpecific(*B,MAT_CLASSID,4);   */
1668   {
1669     PetscBool flag = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_operator = PETSC_FALSE;
1670     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit",&flag,PETSC_NULL);CHKERRQ(ierr);
1671     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw",&flag_draw,PETSC_NULL);CHKERRQ(ierr);
1672     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_explicit_draw_contour",&flag_contour,PETSC_NULL);CHKERRQ(ierr);
1673     ierr  = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_operator",&flag_operator,PETSC_NULL);CHKERRQ(ierr);
1674     if (flag || flag_draw || flag_contour) {
1675       Mat Bexp_mine = PETSC_NULL,Bexp,FDexp;
1676       MatStructure mstruct;
1677       PetscViewer vdraw,vstdout;
1678       PetscBool flg;
1679       if (flag_operator) {
1680         ierr = MatComputeExplicitOperator(*A,&Bexp_mine);CHKERRQ(ierr);
1681         Bexp = Bexp_mine;
1682       } else {
1683         /* See if the preconditioning matrix can be viewed and added directly */
1684         ierr = PetscTypeCompareAny((PetscObject)*B,&flg,MATSEQAIJ,MATMPIAIJ,MATSEQDENSE,MATMPIDENSE,MATSEQBAIJ,MATMPIBAIJ,MATSEQSBAIJ,MATMPIBAIJ,"");CHKERRQ(ierr);
1685         if (flg) Bexp = *B;
1686         else {
1687           /* If the "preconditioning" matrix is itself MATSHELL or some other type without direct support */
1688           ierr = MatComputeExplicitOperator(*B,&Bexp_mine);CHKERRQ(ierr);
1689           Bexp = Bexp_mine;
1690         }
1691       }
1692       ierr = MatConvert(Bexp,MATSAME,MAT_INITIAL_MATRIX,&FDexp);CHKERRQ(ierr);
1693       ierr = SNESDefaultComputeJacobian(snes,X,&FDexp,&FDexp,&mstruct,NULL);CHKERRQ(ierr);
1694       ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&vstdout);CHKERRQ(ierr);
1695       if (flag_draw || flag_contour) {
1696         ierr = PetscViewerDrawOpen(((PetscObject)snes)->comm,0,"Explicit Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr);
1697         if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);}
1698       } else vdraw = PETSC_NULL;
1699       ierr = PetscViewerASCIIPrintf(vstdout,"Explicit %s\n",flag_operator?"Jacobian":"preconditioning Jacobian");CHKERRQ(ierr);
1700       if (flag) {ierr = MatView(Bexp,vstdout);CHKERRQ(ierr);}
1701       if (vdraw) {ierr = MatView(Bexp,vdraw);CHKERRQ(ierr);}
1702       ierr = PetscViewerASCIIPrintf(vstdout,"Finite difference Jacobian\n");CHKERRQ(ierr);
1703       if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);}
1704       if (vdraw) {ierr = MatView(FDexp,vdraw);CHKERRQ(ierr);}
1705       ierr = MatAYPX(FDexp,-1.0,Bexp,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1706       ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian\n");CHKERRQ(ierr);
1707       if (flag) {ierr = MatView(FDexp,vstdout);CHKERRQ(ierr);}
1708       if (vdraw) {              /* Always use contour for the difference */
1709         ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);
1710         ierr = MatView(FDexp,vdraw);CHKERRQ(ierr);
1711         ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);
1712       }
1713       if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);}
1714       ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr);
1715       ierr = MatDestroy(&Bexp_mine);CHKERRQ(ierr);
1716       ierr = MatDestroy(&FDexp);CHKERRQ(ierr);
1717     }
1718   }
1719   {
1720     PetscBool flag = PETSC_FALSE,flag_display = PETSC_FALSE,flag_draw = PETSC_FALSE,flag_contour = PETSC_FALSE,flag_threshold = PETSC_FALSE;
1721     PetscReal threshold_atol = PETSC_SQRT_MACHINE_EPSILON,threshold_rtol = 10*PETSC_SQRT_MACHINE_EPSILON;
1722     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring",&flag,PETSC_NULL);CHKERRQ(ierr);
1723     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_display",&flag_display,PETSC_NULL);CHKERRQ(ierr);
1724     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw",&flag_draw,PETSC_NULL);CHKERRQ(ierr);
1725     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_draw_contour",&flag_contour,PETSC_NULL);CHKERRQ(ierr);
1726     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold",&flag_threshold,PETSC_NULL);CHKERRQ(ierr);
1727     ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_rtol",&threshold_rtol,PETSC_NULL);CHKERRQ(ierr);
1728     ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_compare_coloring_threshold_atol",&threshold_atol,PETSC_NULL);CHKERRQ(ierr);
1729     if (flag || flag_display || flag_draw || flag_contour || flag_threshold) {
1730       Mat Bfd;
1731       MatStructure mstruct;
1732       PetscViewer vdraw,vstdout;
1733       ISColoring iscoloring;
1734       MatFDColoring matfdcoloring;
1735       PetscErrorCode (*func)(SNES,Vec,Vec,void*);
1736       void *funcctx;
1737       PetscReal norm1,norm2,normmax;
1738 
1739       ierr = MatDuplicate(*B,MAT_DO_NOT_COPY_VALUES,&Bfd);CHKERRQ(ierr);
1740       ierr = MatGetColoring(Bfd,MATCOLORINGSL,&iscoloring);CHKERRQ(ierr);
1741       ierr = MatFDColoringCreate(Bfd,iscoloring,&matfdcoloring);CHKERRQ(ierr);
1742       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
1743 
1744       /* This method of getting the function is currently unreliable since it doesn't work for DM local functions. */
1745       ierr = SNESGetFunction(snes,PETSC_NULL,&func,&funcctx);CHKERRQ(ierr);
1746       ierr = MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode(*)(void))func,funcctx);CHKERRQ(ierr);
1747       ierr = PetscObjectSetOptionsPrefix((PetscObject)matfdcoloring,((PetscObject)snes)->prefix);CHKERRQ(ierr);
1748       ierr = PetscObjectAppendOptionsPrefix((PetscObject)matfdcoloring,"coloring_");CHKERRQ(ierr);
1749       ierr = MatFDColoringSetFromOptions(matfdcoloring);CHKERRQ(ierr);
1750       ierr = MatFDColoringApply(Bfd,matfdcoloring,X,&mstruct,snes);CHKERRQ(ierr);
1751       ierr = MatFDColoringDestroy(&matfdcoloring);CHKERRQ(ierr);
1752 
1753       ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&vstdout);CHKERRQ(ierr);
1754       if (flag_draw || flag_contour) {
1755         ierr = PetscViewerDrawOpen(((PetscObject)snes)->comm,0,"Colored Jacobians",PETSC_DECIDE,PETSC_DECIDE,300,300,&vdraw);CHKERRQ(ierr);
1756         if (flag_contour) {ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);}
1757       } else vdraw = PETSC_NULL;
1758       ierr = PetscViewerASCIIPrintf(vstdout,"Explicit preconditioning Jacobian\n");CHKERRQ(ierr);
1759       if (flag_display) {ierr = MatView(*B,vstdout);CHKERRQ(ierr);}
1760       if (vdraw) {ierr = MatView(*B,vdraw);CHKERRQ(ierr);}
1761       ierr = PetscViewerASCIIPrintf(vstdout,"Colored Finite difference Jacobian\n");CHKERRQ(ierr);
1762       if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);}
1763       if (vdraw) {ierr = MatView(Bfd,vdraw);CHKERRQ(ierr);}
1764       ierr = MatAYPX(Bfd,-1.0,*B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1765       ierr = MatNorm(Bfd,NORM_1,&norm1);CHKERRQ(ierr);
1766       ierr = MatNorm(Bfd,NORM_FROBENIUS,&norm2);CHKERRQ(ierr);
1767       ierr = MatNorm(Bfd,NORM_MAX,&normmax);CHKERRQ(ierr);
1768       ierr = PetscViewerASCIIPrintf(vstdout,"User-provided matrix minus finite difference Jacobian, norm1=%G normFrob=%G normmax=%G\n",norm1,norm2,normmax);CHKERRQ(ierr);
1769       if (flag_display) {ierr = MatView(Bfd,vstdout);CHKERRQ(ierr);}
1770       if (vdraw) {              /* Always use contour for the difference */
1771         ierr = PetscViewerPushFormat(vdraw,PETSC_VIEWER_DRAW_CONTOUR);CHKERRQ(ierr);
1772         ierr = MatView(Bfd,vdraw);CHKERRQ(ierr);
1773         ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);
1774       }
1775       if (flag_contour) {ierr = PetscViewerPopFormat(vdraw);CHKERRQ(ierr);}
1776 
1777       if (flag_threshold) {
1778         PetscInt bs,rstart,rend,i;
1779         ierr = MatGetBlockSize(*B,&bs);CHKERRQ(ierr);
1780         ierr = MatGetOwnershipRange(*B,&rstart,&rend);CHKERRQ(ierr);
1781         for (i=rstart; i<rend; i++) {
1782           const PetscScalar *ba,*ca;
1783           const PetscInt *bj,*cj;
1784           PetscInt bn,cn,j,maxentrycol = -1,maxdiffcol = -1,maxrdiffcol = -1;
1785           PetscReal maxentry = 0,maxdiff = 0,maxrdiff = 0;
1786           ierr = MatGetRow(*B,i,&bn,&bj,&ba);CHKERRQ(ierr);
1787           ierr = MatGetRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr);
1788           if (bn != cn) SETERRQ(((PetscObject)*A)->comm,PETSC_ERR_PLIB,"Unexpected different nonzero pattern in -snes_compare_coloring_threshold");
1789           for (j=0; j<bn; j++) {
1790             PetscReal rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
1791             if (PetscAbsScalar(ba[j]) > PetscAbs(maxentry)) {
1792               maxentrycol = bj[j];
1793               maxentry = PetscRealPart(ba[j]);
1794             }
1795             if (PetscAbsScalar(ca[j]) > PetscAbs(maxdiff)) {
1796               maxdiffcol = bj[j];
1797               maxdiff = PetscRealPart(ca[j]);
1798             }
1799             if (rdiff > maxrdiff) {
1800               maxrdiffcol = bj[j];
1801               maxrdiff = rdiff;
1802             }
1803           }
1804           if (maxrdiff > 1) {
1805             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);
1806             for (j=0; j<bn; j++) {
1807               PetscReal rdiff;
1808               rdiff = PetscAbsScalar(ca[j]) / (threshold_atol + threshold_rtol*PetscAbsScalar(ba[j]));
1809               if (rdiff > 1) {
1810                 ierr = PetscViewerASCIIPrintf(vstdout," (%D,%G:%G)",bj[j],PetscRealPart(ba[j]),PetscRealPart(ca[j]));CHKERRQ(ierr);
1811               }
1812             }
1813             ierr = PetscViewerASCIIPrintf(vstdout,"\n",i,maxentry,maxdiff,maxrdiff);CHKERRQ(ierr);
1814           }
1815           ierr = MatRestoreRow(*B,i,&bn,&bj,&ba);CHKERRQ(ierr);
1816           ierr = MatRestoreRow(Bfd,i,&cn,&cj,&ca);CHKERRQ(ierr);
1817         }
1818       }
1819       ierr = PetscViewerDestroy(&vdraw);CHKERRQ(ierr);
1820       ierr = MatDestroy(&Bfd);CHKERRQ(ierr);
1821     }
1822   }
1823   PetscFunctionReturn(0);
1824 }
1825 
1826 #undef __FUNCT__
1827 #define __FUNCT__ "SNESSetJacobian"
1828 /*@C
1829    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1830    location to store the matrix.
1831 
1832    Logically Collective on SNES and Mat
1833 
1834    Input Parameters:
1835 +  snes - the SNES context
1836 .  A - Jacobian matrix
1837 .  B - preconditioner matrix (usually same as the Jacobian)
1838 .  func - Jacobian evaluation routine (if PETSC_NULL then SNES retains any previously set value)
1839 -  ctx - [optional] user-defined context for private data for the
1840          Jacobian evaluation routine (may be PETSC_NULL) (if PETSC_NULL then SNES retains any previously set value)
1841 
1842    Calling sequence of func:
1843 $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
1844 
1845 +  x - input vector
1846 .  A - Jacobian matrix
1847 .  B - preconditioner matrix, usually the same as A
1848 .  flag - flag indicating information about the preconditioner matrix
1849    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1850 -  ctx - [optional] user-defined Jacobian context
1851 
1852    Notes:
1853    See KSPSetOperators() for important information about setting the flag
1854    output parameter in the routine func().  Be sure to read this information!
1855 
1856    The routine func() takes Mat * as the matrix arguments rather than Mat.
1857    This allows the Jacobian evaluation routine to replace A and/or B with a
1858    completely new new matrix structure (not just different matrix elements)
1859    when appropriate, for instance, if the nonzero structure is changing
1860    throughout the global iterations.
1861 
1862    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
1863    each matrix.
1864 
1865    If using SNESDefaultComputeJacobianColor() to assemble a Jacobian, the ctx argument
1866    must be a MatFDColoring.
1867 
1868    Other defect-correction schemes can be used by computing a different matrix in place of the Jacobian.  One common
1869    example is to use the "Picard linearization" which only differentiates through the highest order parts of each term.
1870 
1871    Level: beginner
1872 
1873 .keywords: SNES, nonlinear, set, Jacobian, matrix
1874 
1875 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
1876 @*/
1877 PetscErrorCode  SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
1878 {
1879   PetscErrorCode ierr;
1880 
1881   PetscFunctionBegin;
1882   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1883   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
1884   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
1885   if (A) PetscCheckSameComm(snes,1,A,2);
1886   if (B) PetscCheckSameComm(snes,1,B,3);
1887   if (func) snes->ops->computejacobian = func;
1888   if (ctx)  snes->jacP                 = ctx;
1889   if (A) {
1890     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
1891     ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr);
1892     snes->jacobian = A;
1893   }
1894   if (B) {
1895     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
1896     ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr);
1897     snes->jacobian_pre = B;
1898   }
1899   PetscFunctionReturn(0);
1900 }
1901 
1902 #undef __FUNCT__
1903 #define __FUNCT__ "SNESGetJacobian"
1904 /*@C
1905    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1906    provided context for evaluating the Jacobian.
1907 
1908    Not Collective, but Mat object will be parallel if SNES object is
1909 
1910    Input Parameter:
1911 .  snes - the nonlinear solver context
1912 
1913    Output Parameters:
1914 +  A - location to stash Jacobian matrix (or PETSC_NULL)
1915 .  B - location to stash preconditioner matrix (or PETSC_NULL)
1916 .  func - location to put Jacobian function (or PETSC_NULL)
1917 -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1918 
1919    Level: advanced
1920 
1921 .seealso: SNESSetJacobian(), SNESComputeJacobian()
1922 @*/
1923 PetscErrorCode  SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1924 {
1925   PetscFunctionBegin;
1926   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1927   if (A)    *A    = snes->jacobian;
1928   if (B)    *B    = snes->jacobian_pre;
1929   if (func) *func = snes->ops->computejacobian;
1930   if (ctx)  *ctx  = snes->jacP;
1931   PetscFunctionReturn(0);
1932 }
1933 
1934 /* ----- Routines to initialize and destroy a nonlinear solver ---- */
1935 
1936 #undef __FUNCT__
1937 #define __FUNCT__ "SNESSetUp"
1938 /*@
1939    SNESSetUp - Sets up the internal data structures for the later use
1940    of a nonlinear solver.
1941 
1942    Collective on SNES
1943 
1944    Input Parameters:
1945 .  snes - the SNES context
1946 
1947    Notes:
1948    For basic use of the SNES solvers the user need not explicitly call
1949    SNESSetUp(), since these actions will automatically occur during
1950    the call to SNESSolve().  However, if one wishes to control this
1951    phase separately, SNESSetUp() should be called after SNESCreate()
1952    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1953 
1954    Level: advanced
1955 
1956 .keywords: SNES, nonlinear, setup
1957 
1958 .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
1959 @*/
1960 PetscErrorCode  SNESSetUp(SNES snes)
1961 {
1962   PetscErrorCode ierr;
1963 
1964   PetscFunctionBegin;
1965   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1966   if (snes->setupcalled) PetscFunctionReturn(0);
1967 
1968   if (!((PetscObject)snes)->type_name) {
1969     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
1970   }
1971 
1972   if (!snes->vec_func) {
1973     if (snes->vec_rhs) {
1974       ierr = VecDuplicate(snes->vec_rhs,&snes->vec_func);CHKERRQ(ierr);
1975     } else if (snes->vec_sol) {
1976       ierr = VecDuplicate(snes->vec_sol,&snes->vec_func);CHKERRQ(ierr);
1977     } else if (snes->dm) {
1978       ierr = DMCreateGlobalVector(snes->dm,&snes->vec_func);CHKERRQ(ierr);
1979     }
1980   }
1981   if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1982   if (snes->vec_rhs  == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be right hand side vector");
1983 
1984   if (!snes->vec_sol_update /* && snes->vec_sol */) {
1985     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
1986     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
1987   }
1988 
1989   if (!snes->ops->computejacobian && snes->dm) {
1990     Mat J;
1991     ierr = DMCreateMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr);
1992     ierr = SNESSetJacobian(snes,J,J,SNESDMComputeJacobian,PETSC_NULL);CHKERRQ(ierr);
1993     ierr = MatDestroy(&J);CHKERRQ(ierr);
1994   } else if (!snes->jacobian && snes->ops->computejacobian == MatMFFDComputeJacobian) {
1995     Mat J;
1996     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
1997     ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr);
1998     ierr = MatSetFromOptions(J);CHKERRQ(ierr);
1999     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
2000     ierr = MatDestroy(&J);CHKERRQ(ierr);
2001   } else if (snes->dm && snes->mf_operator && !snes->jacobian_pre && !snes->jacobian) {
2002     Mat J,B;
2003     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
2004     ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr);
2005     ierr = MatSetFromOptions(J);CHKERRQ(ierr);
2006     ierr = DMCreateMatrix(snes->dm,MATAIJ,&B);CHKERRQ(ierr);
2007     ierr = SNESSetJacobian(snes,J,B,SNESDMComputeJacobian,snes->funP);CHKERRQ(ierr);
2008     ierr = MatDestroy(&J);CHKERRQ(ierr);
2009     ierr = MatDestroy(&B);CHKERRQ(ierr);
2010   } else if (snes->dm && !snes->jacobian_pre){
2011     Mat J;
2012     ierr = DMCreateMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr);
2013     ierr = SNESSetJacobian(snes,J,J,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2014     ierr = MatDestroy(&J);CHKERRQ(ierr);
2015   }
2016   if (!snes->ops->computefunction && !snes->dm) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must provide a residual function with SNESSetFunction() or call SNESSetDM()");
2017   if (!snes->vec_func) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must provide a vector when calling SNESSetFunction() or call SNESSetDM()");
2018 
2019   if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);}
2020 
2021   if (snes->ops->usercompute && !snes->user) {
2022     ierr = (*snes->ops->usercompute)(snes,(void**)&snes->user);CHKERRQ(ierr);
2023   }
2024 
2025   if (snes->ops->setup) {
2026     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
2027   }
2028 
2029   snes->setupcalled = PETSC_TRUE;
2030   PetscFunctionReturn(0);
2031 }
2032 
2033 #undef __FUNCT__
2034 #define __FUNCT__ "SNESReset"
2035 /*@
2036    SNESReset - Resets a SNES context to the snessetupcalled = 0 state and removes any allocated Vecs and Mats
2037 
2038    Collective on SNES
2039 
2040    Input Parameter:
2041 .  snes - iterative context obtained from SNESCreate()
2042 
2043    Level: intermediate
2044 
2045    Notes: Also calls the application context destroy routine set with SNESSetComputeApplicationContext()
2046 
2047 .keywords: SNES, destroy
2048 
2049 .seealso: SNESCreate(), SNESSetUp(), SNESSolve()
2050 @*/
2051 PetscErrorCode  SNESReset(SNES snes)
2052 {
2053   PetscErrorCode ierr;
2054 
2055   PetscFunctionBegin;
2056   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2057   if (snes->ops->userdestroy && snes->user) {
2058     ierr       = (*snes->ops->userdestroy)((void**)&snes->user);CHKERRQ(ierr);
2059     snes->user = PETSC_NULL;
2060   }
2061   if (snes->pc) {
2062     ierr = SNESReset(snes->pc);CHKERRQ(ierr);
2063   }
2064 
2065   if (snes->ops->reset) {
2066     ierr = (*snes->ops->reset)(snes);CHKERRQ(ierr);
2067   }
2068   if (snes->ksp) {ierr = KSPReset(snes->ksp);CHKERRQ(ierr);}
2069   ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr);
2070   ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr);
2071   ierr = VecDestroy(&snes->vec_sol_update);CHKERRQ(ierr);
2072   ierr = VecDestroy(&snes->vec_func);CHKERRQ(ierr);
2073   ierr = MatDestroy(&snes->jacobian);CHKERRQ(ierr);
2074   ierr = MatDestroy(&snes->jacobian_pre);CHKERRQ(ierr);
2075   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2076   if (snes->vwork) {ierr = VecDestroyVecs(snes->nvwork,&snes->vwork);CHKERRQ(ierr);}
2077   snes->nwork = snes->nvwork = 0;
2078   snes->setupcalled = PETSC_FALSE;
2079   PetscFunctionReturn(0);
2080 }
2081 
2082 #undef __FUNCT__
2083 #define __FUNCT__ "SNESDestroy"
2084 /*@
2085    SNESDestroy - Destroys the nonlinear solver context that was created
2086    with SNESCreate().
2087 
2088    Collective on SNES
2089 
2090    Input Parameter:
2091 .  snes - the SNES context
2092 
2093    Level: beginner
2094 
2095 .keywords: SNES, nonlinear, destroy
2096 
2097 .seealso: SNESCreate(), SNESSolve()
2098 @*/
2099 PetscErrorCode  SNESDestroy(SNES *snes)
2100 {
2101   PetscErrorCode ierr;
2102 
2103   PetscFunctionBegin;
2104   if (!*snes) PetscFunctionReturn(0);
2105   PetscValidHeaderSpecific((*snes),SNES_CLASSID,1);
2106   if (--((PetscObject)(*snes))->refct > 0) {*snes = 0; PetscFunctionReturn(0);}
2107 
2108   ierr = SNESReset((*snes));CHKERRQ(ierr);
2109   ierr = SNESDestroy(&(*snes)->pc);CHKERRQ(ierr);
2110 
2111   /* if memory was published with AMS then destroy it */
2112   ierr = PetscObjectDepublish((*snes));CHKERRQ(ierr);
2113   if ((*snes)->ops->destroy) {ierr = (*((*snes))->ops->destroy)((*snes));CHKERRQ(ierr);}
2114 
2115   ierr = DMDestroy(&(*snes)->dm);CHKERRQ(ierr);
2116   ierr = KSPDestroy(&(*snes)->ksp);CHKERRQ(ierr);
2117 
2118   ierr = PetscFree((*snes)->kspconvctx);CHKERRQ(ierr);
2119   if ((*snes)->ops->convergeddestroy) {
2120     ierr = (*(*snes)->ops->convergeddestroy)((*snes)->cnvP);CHKERRQ(ierr);
2121   }
2122   if ((*snes)->conv_malloc) {
2123     ierr = PetscFree((*snes)->conv_hist);CHKERRQ(ierr);
2124     ierr = PetscFree((*snes)->conv_hist_its);CHKERRQ(ierr);
2125   }
2126   ierr = PetscViewerDestroy(&(*snes)->ls_monitor);CHKERRQ(ierr);
2127   ierr = SNESMonitorCancel((*snes));CHKERRQ(ierr);
2128   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
2129  PetscFunctionReturn(0);
2130 }
2131 
2132 /* ----------- Routines to set solver parameters ---------- */
2133 
2134 #undef __FUNCT__
2135 #define __FUNCT__ "SNESSetLagPreconditioner"
2136 /*@
2137    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
2138 
2139    Logically Collective on SNES
2140 
2141    Input Parameters:
2142 +  snes - the SNES context
2143 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2144          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2145 
2146    Options Database Keys:
2147 .    -snes_lag_preconditioner <lag>
2148 
2149    Notes:
2150    The default is 1
2151    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2152    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
2153 
2154    Level: intermediate
2155 
2156 .keywords: SNES, nonlinear, set, convergence, tolerances
2157 
2158 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2159 
2160 @*/
2161 PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
2162 {
2163   PetscFunctionBegin;
2164   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2165   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2166   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2167   PetscValidLogicalCollectiveInt(snes,lag,2);
2168   snes->lagpreconditioner = lag;
2169   PetscFunctionReturn(0);
2170 }
2171 
2172 #undef __FUNCT__
2173 #define __FUNCT__ "SNESSetGridSequence"
2174 /*@
2175    SNESSetGridSequence - sets the number of steps of grid sequencing that SNES does
2176 
2177    Logically Collective on SNES
2178 
2179    Input Parameters:
2180 +  snes - the SNES context
2181 -  steps - the number of refinements to do, defaults to 0
2182 
2183    Options Database Keys:
2184 .    -snes_grid_sequence <steps>
2185 
2186    Level: intermediate
2187 
2188    Notes:
2189    Use SNESGetSolution() to extract the fine grid solution after grid sequencing.
2190 
2191 .keywords: SNES, nonlinear, set, convergence, tolerances
2192 
2193 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
2194 
2195 @*/
2196 PetscErrorCode  SNESSetGridSequence(SNES snes,PetscInt steps)
2197 {
2198   PetscFunctionBegin;
2199   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2200   PetscValidLogicalCollectiveInt(snes,steps,2);
2201   snes->gridsequence = steps;
2202   PetscFunctionReturn(0);
2203 }
2204 
2205 #undef __FUNCT__
2206 #define __FUNCT__ "SNESGetLagPreconditioner"
2207 /*@
2208    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
2209 
2210    Not Collective
2211 
2212    Input Parameter:
2213 .  snes - the SNES context
2214 
2215    Output Parameter:
2216 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2217          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
2218 
2219    Options Database Keys:
2220 .    -snes_lag_preconditioner <lag>
2221 
2222    Notes:
2223    The default is 1
2224    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2225 
2226    Level: intermediate
2227 
2228 .keywords: SNES, nonlinear, set, convergence, tolerances
2229 
2230 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
2231 
2232 @*/
2233 PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
2234 {
2235   PetscFunctionBegin;
2236   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2237   *lag = snes->lagpreconditioner;
2238   PetscFunctionReturn(0);
2239 }
2240 
2241 #undef __FUNCT__
2242 #define __FUNCT__ "SNESSetLagJacobian"
2243 /*@
2244    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
2245      often the preconditioner is rebuilt.
2246 
2247    Logically Collective on SNES
2248 
2249    Input Parameters:
2250 +  snes - the SNES context
2251 -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2252          the Jacobian is built etc. -2 means rebuild at next chance but then never again
2253 
2254    Options Database Keys:
2255 .    -snes_lag_jacobian <lag>
2256 
2257    Notes:
2258    The default is 1
2259    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2260    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
2261    at the next Newton step but never again (unless it is reset to another value)
2262 
2263    Level: intermediate
2264 
2265 .keywords: SNES, nonlinear, set, convergence, tolerances
2266 
2267 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
2268 
2269 @*/
2270 PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
2271 {
2272   PetscFunctionBegin;
2273   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2274   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
2275   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
2276   PetscValidLogicalCollectiveInt(snes,lag,2);
2277   snes->lagjacobian = lag;
2278   PetscFunctionReturn(0);
2279 }
2280 
2281 #undef __FUNCT__
2282 #define __FUNCT__ "SNESGetLagJacobian"
2283 /*@
2284    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
2285 
2286    Not Collective
2287 
2288    Input Parameter:
2289 .  snes - the SNES context
2290 
2291    Output Parameter:
2292 .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
2293          the Jacobian is built etc.
2294 
2295    Options Database Keys:
2296 .    -snes_lag_jacobian <lag>
2297 
2298    Notes:
2299    The default is 1
2300    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
2301 
2302    Level: intermediate
2303 
2304 .keywords: SNES, nonlinear, set, convergence, tolerances
2305 
2306 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
2307 
2308 @*/
2309 PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
2310 {
2311   PetscFunctionBegin;
2312   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2313   *lag = snes->lagjacobian;
2314   PetscFunctionReturn(0);
2315 }
2316 
2317 #undef __FUNCT__
2318 #define __FUNCT__ "SNESSetTolerances"
2319 /*@
2320    SNESSetTolerances - Sets various parameters used in convergence tests.
2321 
2322    Logically Collective on SNES
2323 
2324    Input Parameters:
2325 +  snes - the SNES context
2326 .  abstol - absolute convergence tolerance
2327 .  rtol - relative convergence tolerance
2328 .  stol -  convergence tolerance in terms of the norm
2329            of the change in the solution between steps
2330 .  maxit - maximum number of iterations
2331 -  maxf - maximum number of function evaluations
2332 
2333    Options Database Keys:
2334 +    -snes_atol <abstol> - Sets abstol
2335 .    -snes_rtol <rtol> - Sets rtol
2336 .    -snes_stol <stol> - Sets stol
2337 .    -snes_max_it <maxit> - Sets maxit
2338 -    -snes_max_funcs <maxf> - Sets maxf
2339 
2340    Notes:
2341    The default maximum number of iterations is 50.
2342    The default maximum number of function evaluations is 1000.
2343 
2344    Level: intermediate
2345 
2346 .keywords: SNES, nonlinear, set, convergence, tolerances
2347 
2348 .seealso: SNESSetTrustRegionTolerance()
2349 @*/
2350 PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
2351 {
2352   PetscFunctionBegin;
2353   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2354   PetscValidLogicalCollectiveReal(snes,abstol,2);
2355   PetscValidLogicalCollectiveReal(snes,rtol,3);
2356   PetscValidLogicalCollectiveReal(snes,stol,4);
2357   PetscValidLogicalCollectiveInt(snes,maxit,5);
2358   PetscValidLogicalCollectiveInt(snes,maxf,6);
2359 
2360   if (abstol != PETSC_DEFAULT) {
2361     if (abstol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",abstol);
2362     snes->abstol = abstol;
2363   }
2364   if (rtol != PETSC_DEFAULT) {
2365     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);
2366     snes->rtol = rtol;
2367   }
2368   if (stol != PETSC_DEFAULT) {
2369     if (stol < 0.0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %G must be non-negative",stol);
2370     snes->xtol = stol;
2371   }
2372   if (maxit != PETSC_DEFAULT) {
2373     if (maxit < 0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
2374     snes->max_its = maxit;
2375   }
2376   if (maxf != PETSC_DEFAULT) {
2377     if (maxf < 0) SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of function evaluations %D must be non-negative",maxf);
2378     snes->max_funcs = maxf;
2379   }
2380   PetscFunctionReturn(0);
2381 }
2382 
2383 #undef __FUNCT__
2384 #define __FUNCT__ "SNESGetTolerances"
2385 /*@
2386    SNESGetTolerances - Gets various parameters used in convergence tests.
2387 
2388    Not Collective
2389 
2390    Input Parameters:
2391 +  snes - the SNES context
2392 .  atol - absolute convergence tolerance
2393 .  rtol - relative convergence tolerance
2394 .  stol -  convergence tolerance in terms of the norm
2395            of the change in the solution between steps
2396 .  maxit - maximum number of iterations
2397 -  maxf - maximum number of function evaluations
2398 
2399    Notes:
2400    The user can specify PETSC_NULL for any parameter that is not needed.
2401 
2402    Level: intermediate
2403 
2404 .keywords: SNES, nonlinear, get, convergence, tolerances
2405 
2406 .seealso: SNESSetTolerances()
2407 @*/
2408 PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
2409 {
2410   PetscFunctionBegin;
2411   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2412   if (atol)  *atol  = snes->abstol;
2413   if (rtol)  *rtol  = snes->rtol;
2414   if (stol)  *stol  = snes->xtol;
2415   if (maxit) *maxit = snes->max_its;
2416   if (maxf)  *maxf  = snes->max_funcs;
2417   PetscFunctionReturn(0);
2418 }
2419 
2420 #undef __FUNCT__
2421 #define __FUNCT__ "SNESSetTrustRegionTolerance"
2422 /*@
2423    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
2424 
2425    Logically Collective on SNES
2426 
2427    Input Parameters:
2428 +  snes - the SNES context
2429 -  tol - tolerance
2430 
2431    Options Database Key:
2432 .  -snes_trtol <tol> - Sets tol
2433 
2434    Level: intermediate
2435 
2436 .keywords: SNES, nonlinear, set, trust region, tolerance
2437 
2438 .seealso: SNESSetTolerances()
2439 @*/
2440 PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
2441 {
2442   PetscFunctionBegin;
2443   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2444   PetscValidLogicalCollectiveReal(snes,tol,2);
2445   snes->deltatol = tol;
2446   PetscFunctionReturn(0);
2447 }
2448 
2449 /*
2450    Duplicate the lg monitors for SNES from KSP; for some reason with
2451    dynamic libraries things don't work under Sun4 if we just use
2452    macros instead of functions
2453 */
2454 #undef __FUNCT__
2455 #define __FUNCT__ "SNESMonitorLG"
2456 PetscErrorCode  SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
2457 {
2458   PetscErrorCode ierr;
2459 
2460   PetscFunctionBegin;
2461   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2462   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
2463   PetscFunctionReturn(0);
2464 }
2465 
2466 #undef __FUNCT__
2467 #define __FUNCT__ "SNESMonitorLGCreate"
2468 PetscErrorCode  SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2469 {
2470   PetscErrorCode ierr;
2471 
2472   PetscFunctionBegin;
2473   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
2474   PetscFunctionReturn(0);
2475 }
2476 
2477 #undef __FUNCT__
2478 #define __FUNCT__ "SNESMonitorLGDestroy"
2479 PetscErrorCode  SNESMonitorLGDestroy(PetscDrawLG *draw)
2480 {
2481   PetscErrorCode ierr;
2482 
2483   PetscFunctionBegin;
2484   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
2485   PetscFunctionReturn(0);
2486 }
2487 
2488 extern PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
2489 #undef __FUNCT__
2490 #define __FUNCT__ "SNESMonitorLGRange"
2491 PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
2492 {
2493   PetscDrawLG      lg;
2494   PetscErrorCode   ierr;
2495   PetscReal        x,y,per;
2496   PetscViewer      v = (PetscViewer)monctx;
2497   static PetscReal prev; /* should be in the context */
2498   PetscDraw        draw;
2499   PetscFunctionBegin;
2500   if (!monctx) {
2501     MPI_Comm    comm;
2502 
2503     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
2504     v      = PETSC_VIEWER_DRAW_(comm);
2505   }
2506   ierr   = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
2507   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2508   ierr   = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2509   ierr   = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
2510   x = (PetscReal) n;
2511   if (rnorm > 0.0) y = log10(rnorm); else y = -15.0;
2512   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2513   if (n < 20 || !(n % 5)) {
2514     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2515   }
2516 
2517   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
2518   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2519   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2520   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
2521   ierr =  SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr);
2522   x = (PetscReal) n;
2523   y = 100.0*per;
2524   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2525   if (n < 20 || !(n % 5)) {
2526     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2527   }
2528 
2529   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
2530   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2531   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2532   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
2533   x = (PetscReal) n;
2534   y = (prev - rnorm)/prev;
2535   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2536   if (n < 20 || !(n % 5)) {
2537     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2538   }
2539 
2540   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
2541   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
2542   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
2543   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
2544   x = (PetscReal) n;
2545   y = (prev - rnorm)/(prev*per);
2546   if (n > 2) { /*skip initial crazy value */
2547     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
2548   }
2549   if (n < 20 || !(n % 5)) {
2550     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
2551   }
2552   prev = rnorm;
2553   PetscFunctionReturn(0);
2554 }
2555 
2556 #undef __FUNCT__
2557 #define __FUNCT__ "SNESMonitorLGRangeCreate"
2558 PetscErrorCode  SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
2559 {
2560   PetscErrorCode ierr;
2561 
2562   PetscFunctionBegin;
2563   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
2564   PetscFunctionReturn(0);
2565 }
2566 
2567 #undef __FUNCT__
2568 #define __FUNCT__ "SNESMonitorLGRangeDestroy"
2569 PetscErrorCode  SNESMonitorLGRangeDestroy(PetscDrawLG *draw)
2570 {
2571   PetscErrorCode ierr;
2572 
2573   PetscFunctionBegin;
2574   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
2575   PetscFunctionReturn(0);
2576 }
2577 
2578 #undef __FUNCT__
2579 #define __FUNCT__ "SNESMonitor"
2580 /*@
2581    SNESMonitor - runs the user provided monitor routines, if they exist
2582 
2583    Collective on SNES
2584 
2585    Input Parameters:
2586 +  snes - nonlinear solver context obtained from SNESCreate()
2587 .  iter - iteration number
2588 -  rnorm - relative norm of the residual
2589 
2590    Notes:
2591    This routine is called by the SNES implementations.
2592    It does not typically need to be called by the user.
2593 
2594    Level: developer
2595 
2596 .seealso: SNESMonitorSet()
2597 @*/
2598 PetscErrorCode  SNESMonitor(SNES snes,PetscInt iter,PetscReal rnorm)
2599 {
2600   PetscErrorCode ierr;
2601   PetscInt       i,n = snes->numbermonitors;
2602 
2603   PetscFunctionBegin;
2604   for (i=0; i<n; i++) {
2605     ierr = (*snes->monitor[i])(snes,iter,rnorm,snes->monitorcontext[i]);CHKERRQ(ierr);
2606   }
2607   PetscFunctionReturn(0);
2608 }
2609 
2610 /* ------------ Routines to set performance monitoring options ----------- */
2611 
2612 #undef __FUNCT__
2613 #define __FUNCT__ "SNESMonitorSet"
2614 /*@C
2615    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
2616    iteration of the nonlinear solver to display the iteration's
2617    progress.
2618 
2619    Logically Collective on SNES
2620 
2621    Input Parameters:
2622 +  snes - the SNES context
2623 .  func - monitoring routine
2624 .  mctx - [optional] user-defined context for private data for the
2625           monitor routine (use PETSC_NULL if no context is desired)
2626 -  monitordestroy - [optional] routine that frees monitor context
2627           (may be PETSC_NULL)
2628 
2629    Calling sequence of func:
2630 $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
2631 
2632 +    snes - the SNES context
2633 .    its - iteration number
2634 .    norm - 2-norm function value (may be estimated)
2635 -    mctx - [optional] monitoring context
2636 
2637    Options Database Keys:
2638 +    -snes_monitor        - sets SNESMonitorDefault()
2639 .    -snes_monitor_draw    - sets line graph monitor,
2640                             uses SNESMonitorLGCreate()
2641 -    -snes_monitor_cancel - cancels all monitors that have
2642                             been hardwired into a code by
2643                             calls to SNESMonitorSet(), but
2644                             does not cancel those set via
2645                             the options database.
2646 
2647    Notes:
2648    Several different monitoring routines may be set by calling
2649    SNESMonitorSet() multiple times; all will be called in the
2650    order in which they were set.
2651 
2652    Fortran notes: Only a single monitor function can be set for each SNES object
2653 
2654    Level: intermediate
2655 
2656 .keywords: SNES, nonlinear, set, monitor
2657 
2658 .seealso: SNESMonitorDefault(), SNESMonitorCancel()
2659 @*/
2660 PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
2661 {
2662   PetscInt       i;
2663   PetscErrorCode ierr;
2664 
2665   PetscFunctionBegin;
2666   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2667   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
2668   for (i=0; i<snes->numbermonitors;i++) {
2669     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) {
2670       if (monitordestroy) {
2671         ierr = (*monitordestroy)(&mctx);CHKERRQ(ierr);
2672       }
2673       PetscFunctionReturn(0);
2674     }
2675   }
2676   snes->monitor[snes->numbermonitors]           = monitor;
2677   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
2678   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
2679   PetscFunctionReturn(0);
2680 }
2681 
2682 #undef __FUNCT__
2683 #define __FUNCT__ "SNESMonitorCancel"
2684 /*@C
2685    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
2686 
2687    Logically Collective on SNES
2688 
2689    Input Parameters:
2690 .  snes - the SNES context
2691 
2692    Options Database Key:
2693 .  -snes_monitor_cancel - cancels all monitors that have been hardwired
2694     into a code by calls to SNESMonitorSet(), but does not cancel those
2695     set via the options database
2696 
2697    Notes:
2698    There is no way to clear one specific monitor from a SNES object.
2699 
2700    Level: intermediate
2701 
2702 .keywords: SNES, nonlinear, set, monitor
2703 
2704 .seealso: SNESMonitorDefault(), SNESMonitorSet()
2705 @*/
2706 PetscErrorCode  SNESMonitorCancel(SNES snes)
2707 {
2708   PetscErrorCode ierr;
2709   PetscInt       i;
2710 
2711   PetscFunctionBegin;
2712   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2713   for (i=0; i<snes->numbermonitors; i++) {
2714     if (snes->monitordestroy[i]) {
2715       ierr = (*snes->monitordestroy[i])(&snes->monitorcontext[i]);CHKERRQ(ierr);
2716     }
2717   }
2718   snes->numbermonitors = 0;
2719   PetscFunctionReturn(0);
2720 }
2721 
2722 #undef __FUNCT__
2723 #define __FUNCT__ "SNESSetConvergenceTest"
2724 /*@C
2725    SNESSetConvergenceTest - Sets the function that is to be used
2726    to test for convergence of the nonlinear iterative solution.
2727 
2728    Logically Collective on SNES
2729 
2730    Input Parameters:
2731 +  snes - the SNES context
2732 .  func - routine to test for convergence
2733 .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
2734 -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
2735 
2736    Calling sequence of func:
2737 $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
2738 
2739 +    snes - the SNES context
2740 .    it - current iteration (0 is the first and is before any Newton step)
2741 .    cctx - [optional] convergence context
2742 .    reason - reason for convergence/divergence
2743 .    xnorm - 2-norm of current iterate
2744 .    gnorm - 2-norm of current step
2745 -    f - 2-norm of function
2746 
2747    Level: advanced
2748 
2749 .keywords: SNES, nonlinear, set, convergence, test
2750 
2751 .seealso: SNESDefaultConverged(), SNESSkipConverged()
2752 @*/
2753 PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2754 {
2755   PetscErrorCode ierr;
2756 
2757   PetscFunctionBegin;
2758   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2759   if (!func) func = SNESSkipConverged;
2760   if (snes->ops->convergeddestroy) {
2761     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
2762   }
2763   snes->ops->converged        = func;
2764   snes->ops->convergeddestroy = destroy;
2765   snes->cnvP                  = cctx;
2766   PetscFunctionReturn(0);
2767 }
2768 
2769 #undef __FUNCT__
2770 #define __FUNCT__ "SNESGetConvergedReason"
2771 /*@
2772    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
2773 
2774    Not Collective
2775 
2776    Input Parameter:
2777 .  snes - the SNES context
2778 
2779    Output Parameter:
2780 .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
2781             manual pages for the individual convergence tests for complete lists
2782 
2783    Level: intermediate
2784 
2785    Notes: Can only be called after the call the SNESSolve() is complete.
2786 
2787 .keywords: SNES, nonlinear, set, convergence, test
2788 
2789 .seealso: SNESSetConvergenceTest(), SNESConvergedReason
2790 @*/
2791 PetscErrorCode  SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
2792 {
2793   PetscFunctionBegin;
2794   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2795   PetscValidPointer(reason,2);
2796   *reason = snes->reason;
2797   PetscFunctionReturn(0);
2798 }
2799 
2800 #undef __FUNCT__
2801 #define __FUNCT__ "SNESSetConvergenceHistory"
2802 /*@
2803    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
2804 
2805    Logically Collective on SNES
2806 
2807    Input Parameters:
2808 +  snes - iterative context obtained from SNESCreate()
2809 .  a   - array to hold history, this array will contain the function norms computed at each step
2810 .  its - integer array holds the number of linear iterations for each solve.
2811 .  na  - size of a and its
2812 -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
2813            else it continues storing new values for new nonlinear solves after the old ones
2814 
2815    Notes:
2816    If 'a' and 'its' are PETSC_NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
2817    default array of length 10000 is allocated.
2818 
2819    This routine is useful, e.g., when running a code for purposes
2820    of accurate performance monitoring, when no I/O should be done
2821    during the section of code that is being timed.
2822 
2823    Level: intermediate
2824 
2825 .keywords: SNES, set, convergence, history
2826 
2827 .seealso: SNESGetConvergenceHistory()
2828 
2829 @*/
2830 PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool  reset)
2831 {
2832   PetscErrorCode ierr;
2833 
2834   PetscFunctionBegin;
2835   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2836   if (na)  PetscValidScalarPointer(a,2);
2837   if (its) PetscValidIntPointer(its,3);
2838   if (na == PETSC_DECIDE || na == PETSC_DEFAULT || !a) {
2839     if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2840     ierr = PetscMalloc(na*sizeof(PetscReal),&a);CHKERRQ(ierr);
2841     ierr = PetscMalloc(na*sizeof(PetscInt),&its);CHKERRQ(ierr);
2842     snes->conv_malloc   = PETSC_TRUE;
2843   }
2844   snes->conv_hist       = a;
2845   snes->conv_hist_its   = its;
2846   snes->conv_hist_max   = na;
2847   snes->conv_hist_len   = 0;
2848   snes->conv_hist_reset = reset;
2849   PetscFunctionReturn(0);
2850 }
2851 
2852 #if defined(PETSC_HAVE_MATLAB_ENGINE)
2853 #include <engine.h>   /* MATLAB include file */
2854 #include <mex.h>      /* MATLAB include file */
2855 EXTERN_C_BEGIN
2856 #undef __FUNCT__
2857 #define __FUNCT__ "SNESGetConvergenceHistoryMatlab"
2858 mxArray *SNESGetConvergenceHistoryMatlab(SNES snes)
2859 {
2860   mxArray        *mat;
2861   PetscInt       i;
2862   PetscReal      *ar;
2863 
2864   PetscFunctionBegin;
2865   mat  = mxCreateDoubleMatrix(snes->conv_hist_len,1,mxREAL);
2866   ar   = (PetscReal*) mxGetData(mat);
2867   for (i=0; i<snes->conv_hist_len; i++) {
2868     ar[i] = snes->conv_hist[i];
2869   }
2870   PetscFunctionReturn(mat);
2871 }
2872 EXTERN_C_END
2873 #endif
2874 
2875 
2876 #undef __FUNCT__
2877 #define __FUNCT__ "SNESGetConvergenceHistory"
2878 /*@C
2879    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
2880 
2881    Not Collective
2882 
2883    Input Parameter:
2884 .  snes - iterative context obtained from SNESCreate()
2885 
2886    Output Parameters:
2887 .  a   - array to hold history
2888 .  its - integer array holds the number of linear iterations (or
2889          negative if not converged) for each solve.
2890 -  na  - size of a and its
2891 
2892    Notes:
2893     The calling sequence for this routine in Fortran is
2894 $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
2895 
2896    This routine is useful, e.g., when running a code for purposes
2897    of accurate performance monitoring, when no I/O should be done
2898    during the section of code that is being timed.
2899 
2900    Level: intermediate
2901 
2902 .keywords: SNES, get, convergence, history
2903 
2904 .seealso: SNESSetConvergencHistory()
2905 
2906 @*/
2907 PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
2908 {
2909   PetscFunctionBegin;
2910   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2911   if (a)   *a   = snes->conv_hist;
2912   if (its) *its = snes->conv_hist_its;
2913   if (na)  *na  = snes->conv_hist_len;
2914   PetscFunctionReturn(0);
2915 }
2916 
2917 #undef __FUNCT__
2918 #define __FUNCT__ "SNESSetUpdate"
2919 /*@C
2920   SNESSetUpdate - Sets the general-purpose update function called
2921   at the beginning of every iteration of the nonlinear solve. Specifically
2922   it is called just before the Jacobian is "evaluated".
2923 
2924   Logically Collective on SNES
2925 
2926   Input Parameters:
2927 . snes - The nonlinear solver context
2928 . func - The function
2929 
2930   Calling sequence of func:
2931 . func (SNES snes, PetscInt step);
2932 
2933 . step - The current step of the iteration
2934 
2935   Level: advanced
2936 
2937   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()
2938         This is not used by most users.
2939 
2940 .keywords: SNES, update
2941 
2942 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
2943 @*/
2944 PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
2945 {
2946   PetscFunctionBegin;
2947   PetscValidHeaderSpecific(snes, SNES_CLASSID,1);
2948   snes->ops->update = func;
2949   PetscFunctionReturn(0);
2950 }
2951 
2952 #undef __FUNCT__
2953 #define __FUNCT__ "SNESDefaultUpdate"
2954 /*@
2955   SNESDefaultUpdate - The default update function which does nothing.
2956 
2957   Not collective
2958 
2959   Input Parameters:
2960 . snes - The nonlinear solver context
2961 . step - The current step of the iteration
2962 
2963   Level: intermediate
2964 
2965 .keywords: SNES, update
2966 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
2967 @*/
2968 PetscErrorCode  SNESDefaultUpdate(SNES snes, PetscInt step)
2969 {
2970   PetscFunctionBegin;
2971   PetscFunctionReturn(0);
2972 }
2973 
2974 #undef __FUNCT__
2975 #define __FUNCT__ "SNESScaleStep_Private"
2976 /*
2977    SNESScaleStep_Private - Scales a step so that its length is less than the
2978    positive parameter delta.
2979 
2980     Input Parameters:
2981 +   snes - the SNES context
2982 .   y - approximate solution of linear system
2983 .   fnorm - 2-norm of current function
2984 -   delta - trust region size
2985 
2986     Output Parameters:
2987 +   gpnorm - predicted function norm at the new point, assuming local
2988     linearization.  The value is zero if the step lies within the trust
2989     region, and exceeds zero otherwise.
2990 -   ynorm - 2-norm of the step
2991 
2992     Note:
2993     For non-trust region methods such as SNESLS, the parameter delta
2994     is set to be the maximum allowable step size.
2995 
2996 .keywords: SNES, nonlinear, scale, step
2997 */
2998 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
2999 {
3000   PetscReal      nrm;
3001   PetscScalar    cnorm;
3002   PetscErrorCode ierr;
3003 
3004   PetscFunctionBegin;
3005   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3006   PetscValidHeaderSpecific(y,VEC_CLASSID,2);
3007   PetscCheckSameComm(snes,1,y,2);
3008 
3009   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
3010   if (nrm > *delta) {
3011      nrm = *delta/nrm;
3012      *gpnorm = (1.0 - nrm)*(*fnorm);
3013      cnorm = nrm;
3014      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
3015      *ynorm = *delta;
3016   } else {
3017      *gpnorm = 0.0;
3018      *ynorm = nrm;
3019   }
3020   PetscFunctionReturn(0);
3021 }
3022 
3023 #undef __FUNCT__
3024 #define __FUNCT__ "SNESSolve"
3025 /*@C
3026    SNESSolve - Solves a nonlinear system F(x) = b.
3027    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
3028 
3029    Collective on SNES
3030 
3031    Input Parameters:
3032 +  snes - the SNES context
3033 .  b - the constant part of the equation F(x) = b, or PETSC_NULL to use zero.
3034 -  x - the solution vector.
3035 
3036    Notes:
3037    The user should initialize the vector,x, with the initial guess
3038    for the nonlinear solve prior to calling SNESSolve.  In particular,
3039    to employ an initial guess of zero, the user should explicitly set
3040    this vector to zero by calling VecSet().
3041 
3042    Level: beginner
3043 
3044 .keywords: SNES, nonlinear, solve
3045 
3046 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetGridSequence(), SNESGetSolution()
3047 @*/
3048 PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
3049 {
3050   PetscErrorCode ierr;
3051   PetscBool      flg;
3052   char           filename[PETSC_MAX_PATH_LEN];
3053   PetscViewer    viewer;
3054   PetscInt       grid;
3055   Vec            xcreated = PETSC_NULL;
3056 
3057   PetscFunctionBegin;
3058   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3059   if (x) PetscValidHeaderSpecific(x,VEC_CLASSID,3);
3060   if (x) PetscCheckSameComm(snes,1,x,3);
3061   if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2);
3062   if (b) PetscCheckSameComm(snes,1,b,2);
3063 
3064   if (!x && snes->dm) {
3065     ierr = DMCreateGlobalVector(snes->dm,&xcreated);CHKERRQ(ierr);
3066     x    = xcreated;
3067   }
3068 
3069   for (grid=0; grid<snes->gridsequence; grid++) {ierr = PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm));CHKERRQ(ierr);}
3070   for (grid=0; grid<snes->gridsequence+1; grid++) {
3071 
3072     /* set solution vector */
3073     if (!grid) {ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);}
3074     ierr = VecDestroy(&snes->vec_sol);CHKERRQ(ierr);
3075     snes->vec_sol = x;
3076     /* set afine vector if provided */
3077     if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
3078     ierr = VecDestroy(&snes->vec_rhs);CHKERRQ(ierr);
3079     snes->vec_rhs = b;
3080 
3081     ierr = SNESSetUp(snes);CHKERRQ(ierr);
3082 
3083     if (!grid) {
3084       if (snes->ops->computeinitialguess) {
3085         ierr = (*snes->ops->computeinitialguess)(snes,snes->vec_sol,snes->initialguessP);CHKERRQ(ierr);
3086       } else if (snes->dm) {
3087         PetscBool ig;
3088         ierr = DMHasInitialGuess(snes->dm,&ig);CHKERRQ(ierr);
3089         if (ig) {
3090           ierr = DMComputeInitialGuess(snes->dm,snes->vec_sol);CHKERRQ(ierr);
3091         }
3092       }
3093     }
3094 
3095     if (snes->conv_hist_reset) snes->conv_hist_len = 0;
3096     snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
3097 
3098     ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
3099     ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
3100     ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
3101     if (snes->domainerror){
3102       snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
3103       snes->domainerror = PETSC_FALSE;
3104     }
3105     if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
3106 
3107     ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
3108     if (flg && !PetscPreLoadingOn) {
3109       ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
3110       ierr = SNESView(snes,viewer);CHKERRQ(ierr);
3111       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3112     }
3113 
3114     flg  = PETSC_FALSE;
3115     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr);
3116     if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
3117     if (snes->printreason) {
3118       ierr = PetscViewerASCIIAddTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),((PetscObject)snes)->tablevel);CHKERRQ(ierr);
3119       if (snes->reason > 0) {
3120         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
3121       } else {
3122         ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
3123       }
3124       ierr = PetscViewerASCIISubtractTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm),((PetscObject)snes)->tablevel);CHKERRQ(ierr);
3125     }
3126 
3127     flg = PETSC_FALSE;
3128     ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view_solution_vtk",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
3129     if (flg) {
3130       PetscViewer viewer;
3131       ierr = PetscViewerCreate(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
3132       ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
3133       ierr = PetscViewerFileSetName(viewer,filename);CHKERRQ(ierr);
3134       ierr = VecView(snes->vec_sol,viewer);CHKERRQ(ierr);
3135       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
3136     }
3137 
3138     if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
3139     if (grid <  snes->gridsequence) {
3140       DM  fine;
3141       Vec xnew;
3142       Mat interp;
3143 
3144       ierr = DMRefine(snes->dm,((PetscObject)snes)->comm,&fine);CHKERRQ(ierr);
3145       ierr = DMCreateInterpolation(snes->dm,fine,&interp,PETSC_NULL);CHKERRQ(ierr);
3146       ierr = DMCreateGlobalVector(fine,&xnew);CHKERRQ(ierr);
3147       ierr = MatInterpolate(interp,x,xnew);CHKERRQ(ierr);
3148       ierr = MatDestroy(&interp);CHKERRQ(ierr);
3149       x    = xnew;
3150 
3151       ierr = SNESReset(snes);CHKERRQ(ierr);
3152       ierr = SNESSetDM(snes,fine);CHKERRQ(ierr);
3153       ierr = DMDestroy(&fine);CHKERRQ(ierr);
3154       ierr = PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm));CHKERRQ(ierr);
3155     }
3156   }
3157   ierr = VecDestroy(&xcreated);CHKERRQ(ierr);
3158   PetscFunctionReturn(0);
3159 }
3160 
3161 /* --------- Internal routines for SNES Package --------- */
3162 
3163 #undef __FUNCT__
3164 #define __FUNCT__ "SNESSetType"
3165 /*@C
3166    SNESSetType - Sets the method for the nonlinear solver.
3167 
3168    Collective on SNES
3169 
3170    Input Parameters:
3171 +  snes - the SNES context
3172 -  type - a known method
3173 
3174    Options Database Key:
3175 .  -snes_type <type> - Sets the method; use -help for a list
3176    of available methods (for instance, ls or tr)
3177 
3178    Notes:
3179    See "petsc/include/petscsnes.h" for available methods (for instance)
3180 +    SNESLS - Newton's method with line search
3181      (systems of nonlinear equations)
3182 .    SNESTR - Newton's method with trust region
3183      (systems of nonlinear equations)
3184 
3185   Normally, it is best to use the SNESSetFromOptions() command and then
3186   set the SNES solver type from the options database rather than by using
3187   this routine.  Using the options database provides the user with
3188   maximum flexibility in evaluating the many nonlinear solvers.
3189   The SNESSetType() routine is provided for those situations where it
3190   is necessary to set the nonlinear solver independently of the command
3191   line or options database.  This might be the case, for example, when
3192   the choice of solver changes during the execution of the program,
3193   and the user's application is taking responsibility for choosing the
3194   appropriate method.
3195 
3196   Level: intermediate
3197 
3198 .keywords: SNES, set, type
3199 
3200 .seealso: SNESType, SNESCreate()
3201 
3202 @*/
3203 PetscErrorCode  SNESSetType(SNES snes,const SNESType type)
3204 {
3205   PetscErrorCode ierr,(*r)(SNES);
3206   PetscBool      match;
3207 
3208   PetscFunctionBegin;
3209   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3210   PetscValidCharPointer(type,2);
3211 
3212   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
3213   if (match) PetscFunctionReturn(0);
3214 
3215   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr);
3216   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
3217   /* Destroy the previous private SNES context */
3218   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
3219   /* Reinitialize function pointers in SNESOps structure */
3220   snes->ops->setup          = 0;
3221   snes->ops->solve          = 0;
3222   snes->ops->view           = 0;
3223   snes->ops->setfromoptions = 0;
3224   snes->ops->destroy        = 0;
3225   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
3226   snes->setupcalled = PETSC_FALSE;
3227   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
3228   ierr = (*r)(snes);CHKERRQ(ierr);
3229 #if defined(PETSC_HAVE_AMS)
3230   if (PetscAMSPublishAll) {
3231     ierr = PetscObjectAMSPublish((PetscObject)snes);CHKERRQ(ierr);
3232   }
3233 #endif
3234   PetscFunctionReturn(0);
3235 }
3236 
3237 
3238 /* --------------------------------------------------------------------- */
3239 #undef __FUNCT__
3240 #define __FUNCT__ "SNESRegisterDestroy"
3241 /*@
3242    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
3243    registered by SNESRegisterDynamic().
3244 
3245    Not Collective
3246 
3247    Level: advanced
3248 
3249 .keywords: SNES, nonlinear, register, destroy
3250 
3251 .seealso: SNESRegisterAll(), SNESRegisterAll()
3252 @*/
3253 PetscErrorCode  SNESRegisterDestroy(void)
3254 {
3255   PetscErrorCode ierr;
3256 
3257   PetscFunctionBegin;
3258   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
3259   SNESRegisterAllCalled = PETSC_FALSE;
3260   PetscFunctionReturn(0);
3261 }
3262 
3263 #undef __FUNCT__
3264 #define __FUNCT__ "SNESGetType"
3265 /*@C
3266    SNESGetType - Gets the SNES method type and name (as a string).
3267 
3268    Not Collective
3269 
3270    Input Parameter:
3271 .  snes - nonlinear solver context
3272 
3273    Output Parameter:
3274 .  type - SNES method (a character string)
3275 
3276    Level: intermediate
3277 
3278 .keywords: SNES, nonlinear, get, type, name
3279 @*/
3280 PetscErrorCode  SNESGetType(SNES snes,const SNESType *type)
3281 {
3282   PetscFunctionBegin;
3283   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3284   PetscValidPointer(type,2);
3285   *type = ((PetscObject)snes)->type_name;
3286   PetscFunctionReturn(0);
3287 }
3288 
3289 #undef __FUNCT__
3290 #define __FUNCT__ "SNESGetSolution"
3291 /*@
3292    SNESGetSolution - Returns the vector where the approximate solution is
3293    stored. This is the fine grid solution when using SNESSetGridSequence().
3294 
3295    Not Collective, but Vec is parallel if SNES is parallel
3296 
3297    Input Parameter:
3298 .  snes - the SNES context
3299 
3300    Output Parameter:
3301 .  x - the solution
3302 
3303    Level: intermediate
3304 
3305 .keywords: SNES, nonlinear, get, solution
3306 
3307 .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
3308 @*/
3309 PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
3310 {
3311   PetscFunctionBegin;
3312   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3313   PetscValidPointer(x,2);
3314   *x = snes->vec_sol;
3315   PetscFunctionReturn(0);
3316 }
3317 
3318 #undef __FUNCT__
3319 #define __FUNCT__ "SNESGetSolutionUpdate"
3320 /*@
3321    SNESGetSolutionUpdate - Returns the vector where the solution update is
3322    stored.
3323 
3324    Not Collective, but Vec is parallel if SNES is parallel
3325 
3326    Input Parameter:
3327 .  snes - the SNES context
3328 
3329    Output Parameter:
3330 .  x - the solution update
3331 
3332    Level: advanced
3333 
3334 .keywords: SNES, nonlinear, get, solution, update
3335 
3336 .seealso: SNESGetSolution(), SNESGetFunction()
3337 @*/
3338 PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
3339 {
3340   PetscFunctionBegin;
3341   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3342   PetscValidPointer(x,2);
3343   *x = snes->vec_sol_update;
3344   PetscFunctionReturn(0);
3345 }
3346 
3347 #undef __FUNCT__
3348 #define __FUNCT__ "SNESGetFunction"
3349 /*@C
3350    SNESGetFunction - Returns the vector where the function is stored.
3351 
3352    Not Collective, but Vec is parallel if SNES is parallel
3353 
3354    Input Parameter:
3355 .  snes - the SNES context
3356 
3357    Output Parameter:
3358 +  r - the function (or PETSC_NULL)
3359 .  func - the function (or PETSC_NULL)
3360 -  ctx - the function context (or PETSC_NULL)
3361 
3362    Level: advanced
3363 
3364 .keywords: SNES, nonlinear, get, function
3365 
3366 .seealso: SNESSetFunction(), SNESGetSolution()
3367 @*/
3368 PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
3369 {
3370   PetscFunctionBegin;
3371   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3372   if (r)    *r    = snes->vec_func;
3373   if (func) *func = snes->ops->computefunction;
3374   if (ctx)  *ctx  = snes->funP;
3375   PetscFunctionReturn(0);
3376 }
3377 
3378 /*@C
3379    SNESGetGS - Returns the GS function and context.
3380 
3381    Input Parameter:
3382 .  snes - the SNES context
3383 
3384    Output Parameter:
3385 +  gsfunc - the function (or PETSC_NULL)
3386 -  ctx    - the function context (or PETSC_NULL)
3387 
3388    Level: advanced
3389 
3390 .keywords: SNES, nonlinear, get, function
3391 
3392 .seealso: SNESSetGS(), SNESGetFunction()
3393 @*/
3394 
3395 #undef __FUNCT__
3396 #define __FUNCT__ "SNESGetGS"
3397 PetscErrorCode SNESGetGS (SNES snes, PetscErrorCode(**func)(SNES, Vec, Vec, void*), void ** ctx)
3398 {
3399   PetscFunctionBegin;
3400   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3401   if (func) *func = snes->ops->computegs;
3402   if (ctx)  *ctx  = snes->funP;
3403   PetscFunctionReturn(0);
3404 }
3405 
3406 #undef __FUNCT__
3407 #define __FUNCT__ "SNESSetOptionsPrefix"
3408 /*@C
3409    SNESSetOptionsPrefix - Sets the prefix used for searching for all
3410    SNES options in the database.
3411 
3412    Logically Collective on SNES
3413 
3414    Input Parameter:
3415 +  snes - the SNES context
3416 -  prefix - the prefix to prepend to all option names
3417 
3418    Notes:
3419    A hyphen (-) must NOT be given at the beginning of the prefix name.
3420    The first character of all runtime options is AUTOMATICALLY the hyphen.
3421 
3422    Level: advanced
3423 
3424 .keywords: SNES, set, options, prefix, database
3425 
3426 .seealso: SNESSetFromOptions()
3427 @*/
3428 PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
3429 {
3430   PetscErrorCode ierr;
3431 
3432   PetscFunctionBegin;
3433   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3434   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
3435   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
3436   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
3437   PetscFunctionReturn(0);
3438 }
3439 
3440 #undef __FUNCT__
3441 #define __FUNCT__ "SNESAppendOptionsPrefix"
3442 /*@C
3443    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
3444    SNES options in the database.
3445 
3446    Logically Collective on SNES
3447 
3448    Input Parameters:
3449 +  snes - the SNES context
3450 -  prefix - the prefix to prepend to all option names
3451 
3452    Notes:
3453    A hyphen (-) must NOT be given at the beginning of the prefix name.
3454    The first character of all runtime options is AUTOMATICALLY the hyphen.
3455 
3456    Level: advanced
3457 
3458 .keywords: SNES, append, options, prefix, database
3459 
3460 .seealso: SNESGetOptionsPrefix()
3461 @*/
3462 PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
3463 {
3464   PetscErrorCode ierr;
3465 
3466   PetscFunctionBegin;
3467   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3468   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
3469   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
3470   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
3471   PetscFunctionReturn(0);
3472 }
3473 
3474 #undef __FUNCT__
3475 #define __FUNCT__ "SNESGetOptionsPrefix"
3476 /*@C
3477    SNESGetOptionsPrefix - Sets the prefix used for searching for all
3478    SNES options in the database.
3479 
3480    Not Collective
3481 
3482    Input Parameter:
3483 .  snes - the SNES context
3484 
3485    Output Parameter:
3486 .  prefix - pointer to the prefix string used
3487 
3488    Notes: On the fortran side, the user should pass in a string 'prefix' of
3489    sufficient length to hold the prefix.
3490 
3491    Level: advanced
3492 
3493 .keywords: SNES, get, options, prefix, database
3494 
3495 .seealso: SNESAppendOptionsPrefix()
3496 @*/
3497 PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
3498 {
3499   PetscErrorCode ierr;
3500 
3501   PetscFunctionBegin;
3502   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3503   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
3504   PetscFunctionReturn(0);
3505 }
3506 
3507 
3508 #undef __FUNCT__
3509 #define __FUNCT__ "SNESRegister"
3510 /*@C
3511   SNESRegister - See SNESRegisterDynamic()
3512 
3513   Level: advanced
3514 @*/
3515 PetscErrorCode  SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
3516 {
3517   char           fullname[PETSC_MAX_PATH_LEN];
3518   PetscErrorCode ierr;
3519 
3520   PetscFunctionBegin;
3521   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
3522   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
3523   PetscFunctionReturn(0);
3524 }
3525 
3526 #undef __FUNCT__
3527 #define __FUNCT__ "SNESTestLocalMin"
3528 PetscErrorCode  SNESTestLocalMin(SNES snes)
3529 {
3530   PetscErrorCode ierr;
3531   PetscInt       N,i,j;
3532   Vec            u,uh,fh;
3533   PetscScalar    value;
3534   PetscReal      norm;
3535 
3536   PetscFunctionBegin;
3537   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
3538   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
3539   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
3540 
3541   /* currently only works for sequential */
3542   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
3543   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
3544   for (i=0; i<N; i++) {
3545     ierr = VecCopy(u,uh);CHKERRQ(ierr);
3546     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
3547     for (j=-10; j<11; j++) {
3548       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
3549       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
3550       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
3551       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
3552       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
3553       value = -value;
3554       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
3555     }
3556   }
3557   ierr = VecDestroy(&uh);CHKERRQ(ierr);
3558   ierr = VecDestroy(&fh);CHKERRQ(ierr);
3559   PetscFunctionReturn(0);
3560 }
3561 
3562 #undef __FUNCT__
3563 #define __FUNCT__ "SNESKSPSetUseEW"
3564 /*@
3565    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
3566    computing relative tolerance for linear solvers within an inexact
3567    Newton method.
3568 
3569    Logically Collective on SNES
3570 
3571    Input Parameters:
3572 +  snes - SNES context
3573 -  flag - PETSC_TRUE or PETSC_FALSE
3574 
3575     Options Database:
3576 +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
3577 .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
3578 .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
3579 .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
3580 .  -snes_ksp_ew_gamma <gamma> - Sets gamma
3581 .  -snes_ksp_ew_alpha <alpha> - Sets alpha
3582 .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
3583 -  -snes_ksp_ew_threshold <threshold> - Sets threshold
3584 
3585    Notes:
3586    Currently, the default is to use a constant relative tolerance for
3587    the inner linear solvers.  Alternatively, one can use the
3588    Eisenstat-Walker method, where the relative convergence tolerance
3589    is reset at each Newton iteration according progress of the nonlinear
3590    solver.
3591 
3592    Level: advanced
3593 
3594    Reference:
3595    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
3596    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
3597 
3598 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
3599 
3600 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
3601 @*/
3602 PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool  flag)
3603 {
3604   PetscFunctionBegin;
3605   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3606   PetscValidLogicalCollectiveBool(snes,flag,2);
3607   snes->ksp_ewconv = flag;
3608   PetscFunctionReturn(0);
3609 }
3610 
3611 #undef __FUNCT__
3612 #define __FUNCT__ "SNESKSPGetUseEW"
3613 /*@
3614    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
3615    for computing relative tolerance for linear solvers within an
3616    inexact Newton method.
3617 
3618    Not Collective
3619 
3620    Input Parameter:
3621 .  snes - SNES context
3622 
3623    Output Parameter:
3624 .  flag - PETSC_TRUE or PETSC_FALSE
3625 
3626    Notes:
3627    Currently, the default is to use a constant relative tolerance for
3628    the inner linear solvers.  Alternatively, one can use the
3629    Eisenstat-Walker method, where the relative convergence tolerance
3630    is reset at each Newton iteration according progress of the nonlinear
3631    solver.
3632 
3633    Level: advanced
3634 
3635    Reference:
3636    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
3637    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
3638 
3639 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
3640 
3641 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
3642 @*/
3643 PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
3644 {
3645   PetscFunctionBegin;
3646   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3647   PetscValidPointer(flag,2);
3648   *flag = snes->ksp_ewconv;
3649   PetscFunctionReturn(0);
3650 }
3651 
3652 #undef __FUNCT__
3653 #define __FUNCT__ "SNESKSPSetParametersEW"
3654 /*@
3655    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
3656    convergence criteria for the linear solvers within an inexact
3657    Newton method.
3658 
3659    Logically Collective on SNES
3660 
3661    Input Parameters:
3662 +    snes - SNES context
3663 .    version - version 1, 2 (default is 2) or 3
3664 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
3665 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
3666 .    gamma - multiplicative factor for version 2 rtol computation
3667              (0 <= gamma2 <= 1)
3668 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
3669 .    alpha2 - power for safeguard
3670 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
3671 
3672    Note:
3673    Version 3 was contributed by Luis Chacon, June 2006.
3674 
3675    Use PETSC_DEFAULT to retain the default for any of the parameters.
3676 
3677    Level: advanced
3678 
3679    Reference:
3680    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
3681    inexact Newton method", Utah State University Math. Stat. Dept. Res.
3682    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
3683 
3684 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
3685 
3686 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
3687 @*/
3688 PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
3689 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
3690 {
3691   SNESKSPEW *kctx;
3692   PetscFunctionBegin;
3693   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3694   kctx = (SNESKSPEW*)snes->kspconvctx;
3695   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
3696   PetscValidLogicalCollectiveInt(snes,version,2);
3697   PetscValidLogicalCollectiveReal(snes,rtol_0,3);
3698   PetscValidLogicalCollectiveReal(snes,rtol_max,4);
3699   PetscValidLogicalCollectiveReal(snes,gamma,5);
3700   PetscValidLogicalCollectiveReal(snes,alpha,6);
3701   PetscValidLogicalCollectiveReal(snes,alpha2,7);
3702   PetscValidLogicalCollectiveReal(snes,threshold,8);
3703 
3704   if (version != PETSC_DEFAULT)   kctx->version   = version;
3705   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
3706   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
3707   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
3708   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
3709   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
3710   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
3711 
3712   if (kctx->version < 1 || kctx->version > 3) {
3713     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
3714   }
3715   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
3716     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
3717   }
3718   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
3719     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
3720   }
3721   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
3722     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
3723   }
3724   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
3725     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
3726   }
3727   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
3728     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
3729   }
3730   PetscFunctionReturn(0);
3731 }
3732 
3733 #undef __FUNCT__
3734 #define __FUNCT__ "SNESKSPGetParametersEW"
3735 /*@
3736    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
3737    convergence criteria for the linear solvers within an inexact
3738    Newton method.
3739 
3740    Not Collective
3741 
3742    Input Parameters:
3743      snes - SNES context
3744 
3745    Output Parameters:
3746 +    version - version 1, 2 (default is 2) or 3
3747 .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
3748 .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
3749 .    gamma - multiplicative factor for version 2 rtol computation
3750              (0 <= gamma2 <= 1)
3751 .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
3752 .    alpha2 - power for safeguard
3753 -    threshold - threshold for imposing safeguard (0 < threshold < 1)
3754 
3755    Level: advanced
3756 
3757 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
3758 
3759 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
3760 @*/
3761 PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
3762 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
3763 {
3764   SNESKSPEW *kctx;
3765   PetscFunctionBegin;
3766   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3767   kctx = (SNESKSPEW*)snes->kspconvctx;
3768   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
3769   if(version)   *version   = kctx->version;
3770   if(rtol_0)    *rtol_0    = kctx->rtol_0;
3771   if(rtol_max)  *rtol_max  = kctx->rtol_max;
3772   if(gamma)     *gamma     = kctx->gamma;
3773   if(alpha)     *alpha     = kctx->alpha;
3774   if(alpha2)    *alpha2    = kctx->alpha2;
3775   if(threshold) *threshold = kctx->threshold;
3776   PetscFunctionReturn(0);
3777 }
3778 
3779 #undef __FUNCT__
3780 #define __FUNCT__ "SNESKSPEW_PreSolve"
3781 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
3782 {
3783   PetscErrorCode ierr;
3784   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
3785   PetscReal      rtol=PETSC_DEFAULT,stol;
3786 
3787   PetscFunctionBegin;
3788   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
3789   if (!snes->iter) { /* first time in, so use the original user rtol */
3790     rtol = kctx->rtol_0;
3791   } else {
3792     if (kctx->version == 1) {
3793       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
3794       if (rtol < 0.0) rtol = -rtol;
3795       stol = pow(kctx->rtol_last,kctx->alpha2);
3796       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
3797     } else if (kctx->version == 2) {
3798       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
3799       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
3800       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
3801     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
3802       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
3803       /* safeguard: avoid sharp decrease of rtol */
3804       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
3805       stol = PetscMax(rtol,stol);
3806       rtol = PetscMin(kctx->rtol_0,stol);
3807       /* safeguard: avoid oversolving */
3808       stol = kctx->gamma*(snes->ttol)/snes->norm;
3809       stol = PetscMax(rtol,stol);
3810       rtol = PetscMin(kctx->rtol_0,stol);
3811     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
3812   }
3813   /* safeguard: avoid rtol greater than one */
3814   rtol = PetscMin(rtol,kctx->rtol_max);
3815   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
3816   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
3817   PetscFunctionReturn(0);
3818 }
3819 
3820 #undef __FUNCT__
3821 #define __FUNCT__ "SNESKSPEW_PostSolve"
3822 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
3823 {
3824   PetscErrorCode ierr;
3825   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
3826   PCSide         pcside;
3827   Vec            lres;
3828 
3829   PetscFunctionBegin;
3830   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
3831   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
3832   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
3833   if (kctx->version == 1) {
3834     ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr);
3835     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
3836       /* KSP residual is true linear residual */
3837       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
3838     } else {
3839       /* KSP residual is preconditioned residual */
3840       /* compute true linear residual norm */
3841       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
3842       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
3843       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
3844       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
3845       ierr = VecDestroy(&lres);CHKERRQ(ierr);
3846     }
3847   }
3848   PetscFunctionReturn(0);
3849 }
3850 
3851 #undef __FUNCT__
3852 #define __FUNCT__ "SNES_KSPSolve"
3853 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
3854 {
3855   PetscErrorCode ierr;
3856 
3857   PetscFunctionBegin;
3858   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
3859   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
3860   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
3861   PetscFunctionReturn(0);
3862 }
3863 
3864 #undef __FUNCT__
3865 #define __FUNCT__ "SNESSetDM"
3866 /*@
3867    SNESSetDM - Sets the DM that may be used by some preconditioners
3868 
3869    Logically Collective on SNES
3870 
3871    Input Parameters:
3872 +  snes - the preconditioner context
3873 -  dm - the dm
3874 
3875    Level: intermediate
3876 
3877 
3878 .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
3879 @*/
3880 PetscErrorCode  SNESSetDM(SNES snes,DM dm)
3881 {
3882   PetscErrorCode ierr;
3883   KSP            ksp;
3884 
3885   PetscFunctionBegin;
3886   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3887   if (dm) {ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);}
3888   ierr = DMDestroy(&snes->dm);CHKERRQ(ierr);
3889   snes->dm = dm;
3890   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
3891   ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr);
3892   ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr);
3893   if (snes->pc) {
3894     ierr = SNESSetDM(snes->pc, snes->dm);CHKERRQ(ierr);
3895   }
3896   PetscFunctionReturn(0);
3897 }
3898 
3899 #undef __FUNCT__
3900 #define __FUNCT__ "SNESGetDM"
3901 /*@
3902    SNESGetDM - Gets the DM that may be used by some preconditioners
3903 
3904    Not Collective but DM obtained is parallel on SNES
3905 
3906    Input Parameter:
3907 . snes - the preconditioner context
3908 
3909    Output Parameter:
3910 .  dm - the dm
3911 
3912    Level: intermediate
3913 
3914 
3915 .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
3916 @*/
3917 PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
3918 {
3919   PetscFunctionBegin;
3920   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3921   *dm = snes->dm;
3922   PetscFunctionReturn(0);
3923 }
3924 
3925 #undef __FUNCT__
3926 #define __FUNCT__ "SNESSetPC"
3927 /*@
3928   SNESSetPC - Sets the nonlinear preconditioner to be used.
3929 
3930   Collective on SNES
3931 
3932   Input Parameters:
3933 + snes - iterative context obtained from SNESCreate()
3934 - pc   - the preconditioner object
3935 
3936   Notes:
3937   Use SNESGetPC() to retrieve the preconditioner context (for example,
3938   to configure it using the API).
3939 
3940   Level: developer
3941 
3942 .keywords: SNES, set, precondition
3943 .seealso: SNESGetPC()
3944 @*/
3945 PetscErrorCode SNESSetPC(SNES snes, SNES pc)
3946 {
3947   PetscErrorCode ierr;
3948 
3949   PetscFunctionBegin;
3950   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
3951   PetscValidHeaderSpecific(pc, SNES_CLASSID, 2);
3952   PetscCheckSameComm(snes, 1, pc, 2);
3953   ierr = PetscObjectReference((PetscObject) pc);CHKERRQ(ierr);
3954   ierr = SNESDestroy(&snes->pc);CHKERRQ(ierr);
3955   snes->pc = pc;
3956   ierr = PetscLogObjectParent(snes, snes->pc);CHKERRQ(ierr);
3957   PetscFunctionReturn(0);
3958 }
3959 
3960 #undef __FUNCT__
3961 #define __FUNCT__ "SNESGetPC"
3962 /*@
3963   SNESGetPC - Returns a pointer to the nonlinear preconditioning context set with SNESSetPC().
3964 
3965   Not Collective
3966 
3967   Input Parameter:
3968 . snes - iterative context obtained from SNESCreate()
3969 
3970   Output Parameter:
3971 . pc - preconditioner context
3972 
3973   Level: developer
3974 
3975 .keywords: SNES, get, preconditioner
3976 .seealso: SNESSetPC()
3977 @*/
3978 PetscErrorCode SNESGetPC(SNES snes, SNES *pc)
3979 {
3980   PetscErrorCode ierr;
3981 
3982   PetscFunctionBegin;
3983   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
3984   PetscValidPointer(pc, 2);
3985   if (!snes->pc) {
3986     ierr = SNESCreate(((PetscObject) snes)->comm, &snes->pc);CHKERRQ(ierr);
3987     ierr = PetscObjectIncrementTabLevel((PetscObject) snes->pc, (PetscObject) snes, 1);CHKERRQ(ierr);
3988     ierr = PetscLogObjectParent(snes, snes->pc);CHKERRQ(ierr);
3989   }
3990   *pc = snes->pc;
3991   PetscFunctionReturn(0);
3992 }
3993 
3994 #if defined(PETSC_HAVE_MATLAB_ENGINE)
3995 #include <mex.h>
3996 
3997 typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
3998 
3999 #undef __FUNCT__
4000 #define __FUNCT__ "SNESComputeFunction_Matlab"
4001 /*
4002    SNESComputeFunction_Matlab - Calls the function that has been set with
4003                          SNESSetFunctionMatlab().
4004 
4005    Collective on SNES
4006 
4007    Input Parameters:
4008 +  snes - the SNES context
4009 -  x - input vector
4010 
4011    Output Parameter:
4012 .  y - function vector, as set by SNESSetFunction()
4013 
4014    Notes:
4015    SNESComputeFunction() is typically used within nonlinear solvers
4016    implementations, so most users would not generally call this routine
4017    themselves.
4018 
4019    Level: developer
4020 
4021 .keywords: SNES, nonlinear, compute, function
4022 
4023 .seealso: SNESSetFunction(), SNESGetFunction()
4024 */
4025 PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
4026 {
4027   PetscErrorCode    ierr;
4028   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
4029   int               nlhs = 1,nrhs = 5;
4030   mxArray	    *plhs[1],*prhs[5];
4031   long long int     lx = 0,ly = 0,ls = 0;
4032 
4033   PetscFunctionBegin;
4034   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4035   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
4036   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
4037   PetscCheckSameComm(snes,1,x,2);
4038   PetscCheckSameComm(snes,1,y,3);
4039 
4040   /* call Matlab function in ctx with arguments x and y */
4041 
4042   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
4043   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
4044   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
4045   prhs[0] =  mxCreateDoubleScalar((double)ls);
4046   prhs[1] =  mxCreateDoubleScalar((double)lx);
4047   prhs[2] =  mxCreateDoubleScalar((double)ly);
4048   prhs[3] =  mxCreateString(sctx->funcname);
4049   prhs[4] =  sctx->ctx;
4050   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");CHKERRQ(ierr);
4051   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4052   mxDestroyArray(prhs[0]);
4053   mxDestroyArray(prhs[1]);
4054   mxDestroyArray(prhs[2]);
4055   mxDestroyArray(prhs[3]);
4056   mxDestroyArray(plhs[0]);
4057   PetscFunctionReturn(0);
4058 }
4059 
4060 
4061 #undef __FUNCT__
4062 #define __FUNCT__ "SNESSetFunctionMatlab"
4063 /*
4064    SNESSetFunctionMatlab - Sets the function evaluation routine and function
4065    vector for use by the SNES routines in solving systems of nonlinear
4066    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
4067 
4068    Logically Collective on SNES
4069 
4070    Input Parameters:
4071 +  snes - the SNES context
4072 .  r - vector to store function value
4073 -  func - function evaluation routine
4074 
4075    Calling sequence of func:
4076 $    func (SNES snes,Vec x,Vec f,void *ctx);
4077 
4078 
4079    Notes:
4080    The Newton-like methods typically solve linear systems of the form
4081 $      f'(x) x = -f(x),
4082    where f'(x) denotes the Jacobian matrix and f(x) is the function.
4083 
4084    Level: beginner
4085 
4086 .keywords: SNES, nonlinear, set, function
4087 
4088 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4089 */
4090 PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *func,mxArray *ctx)
4091 {
4092   PetscErrorCode    ierr;
4093   SNESMatlabContext *sctx;
4094 
4095   PetscFunctionBegin;
4096   /* currently sctx is memory bleed */
4097   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
4098   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4099   /*
4100      This should work, but it doesn't
4101   sctx->ctx = ctx;
4102   mexMakeArrayPersistent(sctx->ctx);
4103   */
4104   sctx->ctx = mxDuplicateArray(ctx);
4105   ierr = SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);CHKERRQ(ierr);
4106   PetscFunctionReturn(0);
4107 }
4108 
4109 #undef __FUNCT__
4110 #define __FUNCT__ "SNESComputeJacobian_Matlab"
4111 /*
4112    SNESComputeJacobian_Matlab - Calls the function that has been set with
4113                          SNESSetJacobianMatlab().
4114 
4115    Collective on SNES
4116 
4117    Input Parameters:
4118 +  snes - the SNES context
4119 .  x - input vector
4120 .  A, B - the matrices
4121 -  ctx - user context
4122 
4123    Output Parameter:
4124 .  flag - structure of the matrix
4125 
4126    Level: developer
4127 
4128 .keywords: SNES, nonlinear, compute, function
4129 
4130 .seealso: SNESSetFunction(), SNESGetFunction()
4131 @*/
4132 PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag, void *ctx)
4133 {
4134   PetscErrorCode    ierr;
4135   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
4136   int               nlhs = 2,nrhs = 6;
4137   mxArray	    *plhs[2],*prhs[6];
4138   long long int     lx = 0,lA = 0,ls = 0, lB = 0;
4139 
4140   PetscFunctionBegin;
4141   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4142   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
4143 
4144   /* call Matlab function in ctx with arguments x and y */
4145 
4146   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
4147   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
4148   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
4149   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
4150   prhs[0] =  mxCreateDoubleScalar((double)ls);
4151   prhs[1] =  mxCreateDoubleScalar((double)lx);
4152   prhs[2] =  mxCreateDoubleScalar((double)lA);
4153   prhs[3] =  mxCreateDoubleScalar((double)lB);
4154   prhs[4] =  mxCreateString(sctx->funcname);
4155   prhs[5] =  sctx->ctx;
4156   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");CHKERRQ(ierr);
4157   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4158   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
4159   mxDestroyArray(prhs[0]);
4160   mxDestroyArray(prhs[1]);
4161   mxDestroyArray(prhs[2]);
4162   mxDestroyArray(prhs[3]);
4163   mxDestroyArray(prhs[4]);
4164   mxDestroyArray(plhs[0]);
4165   mxDestroyArray(plhs[1]);
4166   PetscFunctionReturn(0);
4167 }
4168 
4169 
4170 #undef __FUNCT__
4171 #define __FUNCT__ "SNESSetJacobianMatlab"
4172 /*
4173    SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
4174    vector for use by the SNES routines in solving systems of nonlinear
4175    equations from MATLAB. Here the function is a string containing the name of a MATLAB function
4176 
4177    Logically Collective on SNES
4178 
4179    Input Parameters:
4180 +  snes - the SNES context
4181 .  A,B - Jacobian matrices
4182 .  func - function evaluation routine
4183 -  ctx - user context
4184 
4185    Calling sequence of func:
4186 $    flag = func (SNES snes,Vec x,Mat A,Mat B,void *ctx);
4187 
4188 
4189    Level: developer
4190 
4191 .keywords: SNES, nonlinear, set, function
4192 
4193 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4194 */
4195 PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *func,mxArray *ctx)
4196 {
4197   PetscErrorCode    ierr;
4198   SNESMatlabContext *sctx;
4199 
4200   PetscFunctionBegin;
4201   /* currently sctx is memory bleed */
4202   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
4203   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4204   /*
4205      This should work, but it doesn't
4206   sctx->ctx = ctx;
4207   mexMakeArrayPersistent(sctx->ctx);
4208   */
4209   sctx->ctx = mxDuplicateArray(ctx);
4210   ierr = SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
4211   PetscFunctionReturn(0);
4212 }
4213 
4214 #undef __FUNCT__
4215 #define __FUNCT__ "SNESMonitor_Matlab"
4216 /*
4217    SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
4218 
4219    Collective on SNES
4220 
4221 .seealso: SNESSetFunction(), SNESGetFunction()
4222 @*/
4223 PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
4224 {
4225   PetscErrorCode  ierr;
4226   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
4227   int             nlhs = 1,nrhs = 6;
4228   mxArray	  *plhs[1],*prhs[6];
4229   long long int   lx = 0,ls = 0;
4230   Vec             x=snes->vec_sol;
4231 
4232   PetscFunctionBegin;
4233   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
4234 
4235   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
4236   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
4237   prhs[0] =  mxCreateDoubleScalar((double)ls);
4238   prhs[1] =  mxCreateDoubleScalar((double)it);
4239   prhs[2] =  mxCreateDoubleScalar((double)fnorm);
4240   prhs[3] =  mxCreateDoubleScalar((double)lx);
4241   prhs[4] =  mxCreateString(sctx->funcname);
4242   prhs[5] =  sctx->ctx;
4243   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");CHKERRQ(ierr);
4244   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
4245   mxDestroyArray(prhs[0]);
4246   mxDestroyArray(prhs[1]);
4247   mxDestroyArray(prhs[2]);
4248   mxDestroyArray(prhs[3]);
4249   mxDestroyArray(prhs[4]);
4250   mxDestroyArray(plhs[0]);
4251   PetscFunctionReturn(0);
4252 }
4253 
4254 
4255 #undef __FUNCT__
4256 #define __FUNCT__ "SNESMonitorSetMatlab"
4257 /*
4258    SNESMonitorSetMatlab - Sets the monitor function from MATLAB
4259 
4260    Level: developer
4261 
4262 .keywords: SNES, nonlinear, set, function
4263 
4264 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
4265 */
4266 PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *func,mxArray *ctx)
4267 {
4268   PetscErrorCode    ierr;
4269   SNESMatlabContext *sctx;
4270 
4271   PetscFunctionBegin;
4272   /* currently sctx is memory bleed */
4273   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
4274   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
4275   /*
4276      This should work, but it doesn't
4277   sctx->ctx = ctx;
4278   mexMakeArrayPersistent(sctx->ctx);
4279   */
4280   sctx->ctx = mxDuplicateArray(ctx);
4281   ierr = SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
4282   PetscFunctionReturn(0);
4283 }
4284 
4285 #endif
4286