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