xref: /petsc/src/snes/interface/snes.c (revision 39d508bb00ddff42cd842b40ab2faf2a7a11a8c1)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
29b94acceSBarry Smith 
37c4f633dSBarry Smith #include "private/snesimpl.h"      /*I "petscsnes.h"  I*/
49b94acceSBarry Smith 
5ace3abfcSBarry Smith PetscBool  SNESRegisterAllCalled = PETSC_FALSE;
68ba1e511SMatthew Knepley PetscFList SNESList              = PETSC_NULL;
78ba1e511SMatthew Knepley 
88ba1e511SMatthew Knepley /* Logging support */
97087cfbeSBarry Smith PetscClassId  SNES_CLASSID;
10166c7f25SBarry Smith PetscLogEvent  SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval;
11a09944afSBarry Smith 
12a09944afSBarry Smith #undef __FUNCT__
13e113a28aSBarry Smith #define __FUNCT__ "SNESSetErrorIfNotConverged"
14e113a28aSBarry Smith /*@
15e113a28aSBarry Smith    SNESSetErrorIfNotConverged - Causes SNESSolve() to generate an error if the solver has not converged.
16e113a28aSBarry Smith 
173f9fe445SBarry Smith    Logically Collective on SNES
18e113a28aSBarry Smith 
19e113a28aSBarry Smith    Input Parameters:
20e113a28aSBarry Smith +  snes - iterative context obtained from SNESCreate()
21e113a28aSBarry Smith -  flg - PETSC_TRUE indicates you want the error generated
22e113a28aSBarry Smith 
23e113a28aSBarry Smith    Options database keys:
24e113a28aSBarry Smith .  -snes_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)
25e113a28aSBarry Smith 
26e113a28aSBarry Smith    Level: intermediate
27e113a28aSBarry Smith 
28e113a28aSBarry Smith    Notes:
29e113a28aSBarry Smith     Normally PETSc continues if a linear solver fails to converge, you can call SNESGetConvergedReason() after a SNESSolve()
30e113a28aSBarry Smith     to determine if it has converged.
31e113a28aSBarry Smith 
32e113a28aSBarry Smith .keywords: SNES, set, initial guess, nonzero
33e113a28aSBarry Smith 
34e113a28aSBarry Smith .seealso: SNESGetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
35e113a28aSBarry Smith @*/
367087cfbeSBarry Smith PetscErrorCode  SNESSetErrorIfNotConverged(SNES snes,PetscBool  flg)
37e113a28aSBarry Smith {
38e113a28aSBarry Smith   PetscFunctionBegin;
39e113a28aSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
40acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(snes,flg,2);
41e113a28aSBarry Smith   snes->errorifnotconverged = flg;
42e113a28aSBarry Smith   PetscFunctionReturn(0);
43e113a28aSBarry Smith }
44e113a28aSBarry Smith 
45e113a28aSBarry Smith #undef __FUNCT__
46e113a28aSBarry Smith #define __FUNCT__ "SNESGetErrorIfNotConverged"
47e113a28aSBarry Smith /*@
48e113a28aSBarry Smith    SNESGetErrorIfNotConverged - Will SNESSolve() generate an error if the solver does not converge?
49e113a28aSBarry Smith 
50e113a28aSBarry Smith    Not Collective
51e113a28aSBarry Smith 
52e113a28aSBarry Smith    Input Parameter:
53e113a28aSBarry Smith .  snes - iterative context obtained from SNESCreate()
54e113a28aSBarry Smith 
55e113a28aSBarry Smith    Output Parameter:
56e113a28aSBarry Smith .  flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE
57e113a28aSBarry Smith 
58e113a28aSBarry Smith    Level: intermediate
59e113a28aSBarry Smith 
60e113a28aSBarry Smith .keywords: SNES, set, initial guess, nonzero
61e113a28aSBarry Smith 
62e113a28aSBarry Smith .seealso:  SNESSetErrorIfNotConverged(), KSPGetErrorIfNotConverged(), KSPSetErrorIFNotConverged()
63e113a28aSBarry Smith @*/
647087cfbeSBarry Smith PetscErrorCode  SNESGetErrorIfNotConverged(SNES snes,PetscBool  *flag)
65e113a28aSBarry Smith {
66e113a28aSBarry Smith   PetscFunctionBegin;
67e113a28aSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
68e113a28aSBarry Smith   PetscValidPointer(flag,2);
69e113a28aSBarry Smith   *flag = snes->errorifnotconverged;
70e113a28aSBarry Smith   PetscFunctionReturn(0);
71e113a28aSBarry Smith }
72e113a28aSBarry Smith 
73e113a28aSBarry Smith #undef __FUNCT__
744936397dSBarry Smith #define __FUNCT__ "SNESSetFunctionDomainError"
75e725d27bSBarry Smith /*@
764936397dSBarry Smith    SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not
774936397dSBarry Smith      in the functions domain. For example, negative pressure.
784936397dSBarry Smith 
793f9fe445SBarry Smith    Logically Collective on SNES
804936397dSBarry Smith 
814936397dSBarry Smith    Input Parameters:
824936397dSBarry Smith .  SNES - the SNES context
834936397dSBarry Smith 
8428529972SSatish Balay    Level: advanced
854936397dSBarry Smith 
864936397dSBarry Smith .keywords: SNES, view
874936397dSBarry Smith 
884936397dSBarry Smith .seealso: SNESCreate(), SNESSetFunction()
894936397dSBarry Smith @*/
907087cfbeSBarry Smith PetscErrorCode  SNESSetFunctionDomainError(SNES snes)
914936397dSBarry Smith {
924936397dSBarry Smith   PetscFunctionBegin;
930700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
944936397dSBarry Smith   snes->domainerror = PETSC_TRUE;
954936397dSBarry Smith   PetscFunctionReturn(0);
964936397dSBarry Smith }
974936397dSBarry Smith 
984936397dSBarry Smith #undef __FUNCT__
994a2ae208SSatish Balay #define __FUNCT__ "SNESView"
1007e2c5f70SBarry Smith /*@C
1019b94acceSBarry Smith    SNESView - Prints the SNES data structure.
1029b94acceSBarry Smith 
1034c49b128SBarry Smith    Collective on SNES
104fee21e36SBarry Smith 
105c7afd0dbSLois Curfman McInnes    Input Parameters:
106c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
107c7afd0dbSLois Curfman McInnes -  viewer - visualization context
108c7afd0dbSLois Curfman McInnes 
1099b94acceSBarry Smith    Options Database Key:
110c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
1119b94acceSBarry Smith 
1129b94acceSBarry Smith    Notes:
1139b94acceSBarry Smith    The available visualization contexts include
114b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
115b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
116c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
117c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
118c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
1199b94acceSBarry Smith 
1203e081fefSLois Curfman McInnes    The user can open an alternative visualization context with
121b0a32e0cSBarry Smith    PetscViewerASCIIOpen() - output to a specified file.
1229b94acceSBarry Smith 
12336851e7fSLois Curfman McInnes    Level: beginner
12436851e7fSLois Curfman McInnes 
1259b94acceSBarry Smith .keywords: SNES, view
1269b94acceSBarry Smith 
127b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
1289b94acceSBarry Smith @*/
1297087cfbeSBarry Smith PetscErrorCode  SNESView(SNES snes,PetscViewer viewer)
1309b94acceSBarry Smith {
131fa9f3622SBarry Smith   SNESKSPEW           *kctx;
132dfbe8321SBarry Smith   PetscErrorCode      ierr;
13394b7f48cSBarry Smith   KSP                 ksp;
134ace3abfcSBarry Smith   PetscBool           iascii,isstring;
1359b94acceSBarry Smith 
1363a40ed3dSBarry Smith   PetscFunctionBegin;
1370700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1383050cee2SBarry Smith   if (!viewer) {
1397adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
1403050cee2SBarry Smith   }
1410700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
142c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,viewer,2);
14374679c65SBarry Smith 
1442692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1452692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
14632077d6dSBarry Smith   if (iascii) {
147317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)snes,viewer,"SNES Object");CHKERRQ(ierr);
148e7788613SBarry Smith     if (snes->ops->view) {
149b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
150e7788613SBarry Smith       ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);
151b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1520ef38995SBarry Smith     }
15377431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
154a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
15570441072SBarry Smith                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
15677431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
15777431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
1589b94acceSBarry Smith     if (snes->ksp_ewconv) {
159fa9f3622SBarry Smith       kctx = (SNESKSPEW *)snes->kspconvctx;
1609b94acceSBarry Smith       if (kctx) {
16177431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
162a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
163a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
1649b94acceSBarry Smith       }
1659b94acceSBarry Smith     }
166eb1f6c34SBarry Smith     if (snes->lagpreconditioner == -1) {
167eb1f6c34SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioned is never rebuilt\n");CHKERRQ(ierr);
168eb1f6c34SBarry Smith     } else if (snes->lagpreconditioner > 1) {
169eb1f6c34SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Preconditioned is rebuilt every %D new Jacobians\n",snes->lagpreconditioner);CHKERRQ(ierr);
170eb1f6c34SBarry Smith     }
171eb1f6c34SBarry Smith     if (snes->lagjacobian == -1) {
172eb1f6c34SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Jacobian is never rebuilt\n");CHKERRQ(ierr);
173eb1f6c34SBarry Smith     } else if (snes->lagjacobian > 1) {
174eb1f6c34SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Jacobian is rebuilt every %D Newton iterations\n",snes->lagjacobian);CHKERRQ(ierr);
175eb1f6c34SBarry Smith     }
1760f5bd95cSBarry Smith   } else if (isstring) {
177317d6ea6SBarry Smith     const char *type;
178454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
179b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
18019bcc07fSBarry Smith   }
18194b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
182b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
18394b7f48cSBarry Smith   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
184b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1853a40ed3dSBarry Smith   PetscFunctionReturn(0);
1869b94acceSBarry Smith }
1879b94acceSBarry Smith 
18876b2cf59SMatthew Knepley /*
18976b2cf59SMatthew Knepley   We retain a list of functions that also take SNES command
19076b2cf59SMatthew Knepley   line options. These are called at the end SNESSetFromOptions()
19176b2cf59SMatthew Knepley */
19276b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5
193a7cc72afSBarry Smith static PetscInt numberofsetfromoptions;
1946849ba73SBarry Smith static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
19576b2cf59SMatthew Knepley 
196e74ef692SMatthew Knepley #undef __FUNCT__
197e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker"
198ac226902SBarry Smith /*@C
19976b2cf59SMatthew Knepley   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
20076b2cf59SMatthew Knepley 
20176b2cf59SMatthew Knepley   Not Collective
20276b2cf59SMatthew Knepley 
20376b2cf59SMatthew Knepley   Input Parameter:
20476b2cf59SMatthew Knepley . snescheck - function that checks for options
20576b2cf59SMatthew Knepley 
20676b2cf59SMatthew Knepley   Level: developer
20776b2cf59SMatthew Knepley 
20876b2cf59SMatthew Knepley .seealso: SNESSetFromOptions()
20976b2cf59SMatthew Knepley @*/
2107087cfbeSBarry Smith PetscErrorCode  SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
21176b2cf59SMatthew Knepley {
21276b2cf59SMatthew Knepley   PetscFunctionBegin;
21376b2cf59SMatthew Knepley   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
214e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
21576b2cf59SMatthew Knepley   }
21676b2cf59SMatthew Knepley   othersetfromoptions[numberofsetfromoptions++] = snescheck;
21776b2cf59SMatthew Knepley   PetscFunctionReturn(0);
21876b2cf59SMatthew Knepley }
21976b2cf59SMatthew Knepley 
2207087cfbeSBarry Smith extern PetscErrorCode  SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
221aa3661deSLisandro Dalcin 
222aa3661deSLisandro Dalcin #undef __FUNCT__
223aa3661deSLisandro Dalcin #define __FUNCT__ "SNESSetUpMatrixFree_Private"
224ace3abfcSBarry Smith static PetscErrorCode SNESSetUpMatrixFree_Private(SNES snes, PetscBool  hasOperator, PetscInt version)
225aa3661deSLisandro Dalcin {
226aa3661deSLisandro Dalcin   Mat            J;
227aa3661deSLisandro Dalcin   KSP            ksp;
228aa3661deSLisandro Dalcin   PC             pc;
229ace3abfcSBarry Smith   PetscBool      match;
230aa3661deSLisandro Dalcin   PetscErrorCode ierr;
231aa3661deSLisandro Dalcin 
232aa3661deSLisandro Dalcin   PetscFunctionBegin;
2330700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
234aa3661deSLisandro Dalcin 
23598613b67SLisandro Dalcin   if(!snes->vec_func && (snes->jacobian || snes->jacobian_pre)) {
23698613b67SLisandro Dalcin     Mat A = snes->jacobian, B = snes->jacobian_pre;
23798613b67SLisandro Dalcin     ierr = MatGetVecs(A ? A : B, PETSC_NULL,&snes->vec_func);CHKERRQ(ierr);
23898613b67SLisandro Dalcin   }
23998613b67SLisandro Dalcin 
240aa3661deSLisandro Dalcin   if (version == 1) {
241aa3661deSLisandro Dalcin     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
24298613b67SLisandro Dalcin     ierr = MatMFFDSetOptionsPrefix(J,((PetscObject)snes)->prefix);CHKERRQ(ierr);
243aa3661deSLisandro Dalcin     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
244aa3661deSLisandro Dalcin   } else if (version == 2) {
245e32f2f54SBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"SNESSetFunction() must be called first");
24665460251SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SCALAR_SINGLE) && !defined(PETSC_USE_SCALAR_MAT_SINGLE) && !defined(PETSC_USE_SCALAR_LONG_DOUBLE) && !defined(PETSC_USE_SCALAR_INT)
247aa3661deSLisandro Dalcin     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_func,&J);CHKERRQ(ierr);
248aa3661deSLisandro Dalcin #else
249e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "matrix-free operator rutines (version 2)");
250aa3661deSLisandro Dalcin #endif
251aa3661deSLisandro Dalcin   } else {
252e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "matrix-free operator rutines, only version 1 and 2");
253aa3661deSLisandro Dalcin   }
254aa3661deSLisandro Dalcin 
255aa3661deSLisandro Dalcin   ierr = PetscInfo1(snes,"Setting default matrix-free operator routines (version %D)\n", version);CHKERRQ(ierr);
256d3462f78SMatthew Knepley   if (hasOperator) {
257aa3661deSLisandro Dalcin     /* This version replaces the user provided Jacobian matrix with a
258aa3661deSLisandro Dalcin        matrix-free version but still employs the user-provided preconditioner matrix. */
259aa3661deSLisandro Dalcin     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
260aa3661deSLisandro Dalcin   } else {
261aa3661deSLisandro Dalcin     /* This version replaces both the user-provided Jacobian and the user-
262aa3661deSLisandro Dalcin        provided preconditioner matrix with the default matrix free version. */
263aa3661deSLisandro Dalcin     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
264aa3661deSLisandro Dalcin     /* Force no preconditioner */
265aa3661deSLisandro Dalcin     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
266aa3661deSLisandro Dalcin     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
267aa3661deSLisandro Dalcin     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&match);CHKERRQ(ierr);
268aa3661deSLisandro Dalcin     if (!match) {
269aa3661deSLisandro Dalcin       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
270aa3661deSLisandro Dalcin       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
271aa3661deSLisandro Dalcin     }
272aa3661deSLisandro Dalcin   }
273aa3661deSLisandro Dalcin   ierr = MatDestroy(J);CHKERRQ(ierr);
274aa3661deSLisandro Dalcin 
275aa3661deSLisandro Dalcin   PetscFunctionReturn(0);
276aa3661deSLisandro Dalcin }
277aa3661deSLisandro Dalcin 
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions"
2809b94acceSBarry Smith /*@
28194b7f48cSBarry Smith    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
2829b94acceSBarry Smith 
283c7afd0dbSLois Curfman McInnes    Collective on SNES
284c7afd0dbSLois Curfman McInnes 
2859b94acceSBarry Smith    Input Parameter:
2869b94acceSBarry Smith .  snes - the SNES context
2879b94acceSBarry Smith 
28836851e7fSLois Curfman McInnes    Options Database Keys:
2896831982aSBarry Smith +  -snes_type <type> - ls, tr, umls, umtr, test
29082738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
29182738288SBarry Smith                 of the change in the solution between steps
29270441072SBarry Smith .  -snes_atol <abstol> - absolute tolerance of residual norm
293b39c3a46SLois Curfman McInnes .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
294b39c3a46SLois Curfman McInnes .  -snes_max_it <max_it> - maximum number of iterations
295b39c3a46SLois Curfman McInnes .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
29650ffb88aSMatthew Knepley .  -snes_max_fail <max_fail> - maximum number of failures
297ddf469c8SBarry Smith .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
298a8054027SBarry Smith .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
299e35cf81dSBarry Smith .  -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild)
300b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
3012492ecdbSBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear
30282738288SBarry Smith                                solver; hence iterations will continue until max_it
3031fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
30482738288SBarry Smith                                of convergence test
305e8105e01SRichard Katz .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
306e8105e01SRichard Katz                                        filename given prints to stdout
307a6570f20SBarry Smith .  -snes_monitor_solution - plots solution at each iteration
308a6570f20SBarry Smith .  -snes_monitor_residual - plots residual (not its norm) at each iteration
309a6570f20SBarry Smith .  -snes_monitor_solution_update - plots update to solution at each iteration
310a6570f20SBarry Smith .  -snes_monitor_draw - plots residual norm at each iteration
311e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
3125968eb51SBarry Smith .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
313fee2055bSBarry Smith -  -snes_converged_reason - print the reason for convergence/divergence after each solve
31482738288SBarry Smith 
31582738288SBarry Smith     Options Database for Eisenstat-Walker method:
316fa9f3622SBarry Smith +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
3174b27c08aSLois Curfman McInnes .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
31836851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
31936851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
32036851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
32136851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
32236851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
32336851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
32482738288SBarry Smith 
32511ca99fdSLois Curfman McInnes    Notes:
32611ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
3270598bfebSBarry Smith    the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A>.
32883e2fdc7SBarry Smith 
32936851e7fSLois Curfman McInnes    Level: beginner
33036851e7fSLois Curfman McInnes 
3319b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
3329b94acceSBarry Smith 
33369ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
3349b94acceSBarry Smith @*/
3357087cfbeSBarry Smith PetscErrorCode  SNESSetFromOptions(SNES snes)
3369b94acceSBarry Smith {
337ace3abfcSBarry Smith   PetscBool               flg,mf,mf_operator;
338aa3661deSLisandro Dalcin   PetscInt                i,indx,lag,mf_version;
339aa3661deSLisandro Dalcin   MatStructure            matflag;
34085385478SLisandro Dalcin   const char              *deft = SNESLS;
34185385478SLisandro Dalcin   const char              *convtests[] = {"default","skip"};
34285385478SLisandro Dalcin   SNESKSPEW               *kctx = NULL;
343e8105e01SRichard Katz   char                    type[256], monfilename[PETSC_MAX_PATH_LEN];
34423d894e5SBarry Smith   PetscViewerASCIIMonitor monviewer;
34585385478SLisandro Dalcin   PetscErrorCode          ierr;
3469b94acceSBarry Smith 
3473a40ed3dSBarry Smith   PetscFunctionBegin;
3480700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
349ca161407SBarry Smith 
350186905e3SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
351cce0b1b2SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
3527adad957SLisandro Dalcin     if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; }
353b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
354d64ed03dSBarry Smith     if (flg) {
355186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
3567adad957SLisandro Dalcin     } else if (!((PetscObject)snes)->type_name) {
357186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
358d64ed03dSBarry Smith     }
35990d69ab7SBarry Smith     /* not used here, but called so will go into help messaage */
360909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
36193c39befSBarry Smith 
36257034d6fSHong Zhang     ierr = PetscOptionsReal("-snes_stol","Stop if step length less than","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
36357034d6fSHong Zhang     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less than","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
364186905e3SBarry Smith 
36557034d6fSHong Zhang     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less than","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
366b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
367b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
36850ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
369ddf469c8SBarry Smith     ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr);
370acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_error_if_not_converged","Generate error if solver does not converge","SNESSetErrorIfNotConverged",snes->errorifnotconverged,&snes->errorifnotconverged,PETSC_NULL);CHKERRQ(ierr);
37185385478SLisandro Dalcin 
372a8054027SBarry Smith     ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr);
373a8054027SBarry Smith     if (flg) {
374a8054027SBarry Smith       ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr);
375a8054027SBarry Smith     }
376e35cf81dSBarry Smith     ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr);
377e35cf81dSBarry Smith     if (flg) {
378e35cf81dSBarry Smith       ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr);
379e35cf81dSBarry Smith     }
380a8054027SBarry Smith 
38185385478SLisandro Dalcin     ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr);
38285385478SLisandro Dalcin     if (flg) {
38385385478SLisandro Dalcin       switch (indx) {
3847f7931b9SBarry Smith       case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break;
3857f7931b9SBarry Smith       case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);    break;
38685385478SLisandro Dalcin       }
38785385478SLisandro Dalcin     }
38885385478SLisandro Dalcin 
389acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",snes->printreason,&snes->printreason,PETSC_NULL);CHKERRQ(ierr);
390186905e3SBarry Smith 
39185385478SLisandro Dalcin     kctx = (SNESKSPEW *)snes->kspconvctx;
39285385478SLisandro Dalcin 
393acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr);
394186905e3SBarry Smith 
395fa9f3622SBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
396fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
397fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
398fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
399fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
400fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
401fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
402186905e3SBarry Smith 
40390d69ab7SBarry Smith     flg  = PETSC_FALSE;
404acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
405a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);}
406eabae89aSBarry Smith 
407a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
408e8105e01SRichard Katz     if (flg) {
409050a712dSBarry Smith       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr);
41023d894e5SBarry Smith       ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
411e8105e01SRichard Katz     }
412eabae89aSBarry Smith 
413b271bb04SBarry Smith     ierr = PetscOptionsString("-snes_monitor_range","Monitor range of elements of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
414b271bb04SBarry Smith     if (flg) {
415b271bb04SBarry Smith       ierr = SNESMonitorSet(snes,SNESMonitorRange,0,0);CHKERRQ(ierr);
416b271bb04SBarry Smith     }
417b271bb04SBarry Smith 
418a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
419eabae89aSBarry Smith     if (flg) {
420050a712dSBarry Smith       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr);
421f1bef1bcSMatthew Knepley       ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr);
422e8105e01SRichard Katz     }
423eabae89aSBarry Smith 
424a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
425eabae89aSBarry Smith     if (flg) {
426050a712dSBarry Smith       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr);
42723d894e5SBarry Smith       ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
428eabae89aSBarry Smith     }
429eabae89aSBarry Smith 
43090d69ab7SBarry Smith     flg  = PETSC_FALSE;
431acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
432a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);}
43390d69ab7SBarry Smith     flg  = PETSC_FALSE;
434acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
435a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);}
43690d69ab7SBarry Smith     flg  = PETSC_FALSE;
437acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
438a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);}
43990d69ab7SBarry Smith     flg  = PETSC_FALSE;
440acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
441a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
44290d69ab7SBarry Smith     flg  = PETSC_FALSE;
443acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_monitor_range_draw","Plot function range at each iteration","SNESMonitorLG",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
444b271bb04SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLGRange,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
445e24b481bSBarry Smith 
44690d69ab7SBarry Smith     flg  = PETSC_FALSE;
447acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
4484b27c08aSLois Curfman McInnes     if (flg) {
449186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
450ae15b995SBarry Smith       ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr);
4519b94acceSBarry Smith     }
452639f9d9dSBarry Smith 
453aa3661deSLisandro Dalcin     mf = mf_operator = PETSC_FALSE;
454aa3661deSLisandro Dalcin     flg = PETSC_FALSE;
455acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_mf_operator","Use a Matrix-Free Jacobian with user-provided preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf_operator,&flg);CHKERRQ(ierr);
456aa3661deSLisandro Dalcin     if (flg && mf_operator) mf = PETSC_TRUE;
457aa3661deSLisandro Dalcin     flg = PETSC_FALSE;
458acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_mf","Use a Matrix-Free Jacobian with no preconditioner matrix","MatCreateSNESMF",PETSC_FALSE,&mf,&flg);CHKERRQ(ierr);
459aa3661deSLisandro Dalcin     if (!flg && mf_operator) mf = PETSC_TRUE;
460aa3661deSLisandro Dalcin     mf_version = 1;
461aa3661deSLisandro Dalcin     ierr = PetscOptionsInt("-snes_mf_version","Matrix-Free routines version 1 or 2","None",mf_version,&mf_version,0);CHKERRQ(ierr);
462aa3661deSLisandro Dalcin 
46376b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
46476b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
46576b2cf59SMatthew Knepley     }
46676b2cf59SMatthew Knepley 
467e7788613SBarry Smith     if (snes->ops->setfromoptions) {
468e7788613SBarry Smith       ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr);
469639f9d9dSBarry Smith     }
4705d973c19SBarry Smith 
4715d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
4725d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)snes);CHKERRQ(ierr);
473b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
4744bbc92c1SBarry Smith 
475aa3661deSLisandro Dalcin   if (mf) { ierr = SNESSetUpMatrixFree_Private(snes, mf_operator, mf_version);CHKERRQ(ierr); }
4761cee3971SBarry Smith 
4771cee3971SBarry Smith   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
478aa3661deSLisandro Dalcin   ierr = KSPGetOperators(snes->ksp,PETSC_NULL,PETSC_NULL,&matflag);
479aa3661deSLisandro Dalcin   ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,matflag);CHKERRQ(ierr);
48085385478SLisandro Dalcin   ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr);
48193993e2dSLois Curfman McInnes 
4823a40ed3dSBarry Smith   PetscFunctionReturn(0);
4839b94acceSBarry Smith }
4849b94acceSBarry Smith 
485a847f771SSatish Balay 
4864a2ae208SSatish Balay #undef __FUNCT__
4874a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
4889b94acceSBarry Smith /*@
4899b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
4909b94acceSBarry Smith    the nonlinear solvers.
4919b94acceSBarry Smith 
4923f9fe445SBarry Smith    Logically Collective on SNES
493fee21e36SBarry Smith 
494c7afd0dbSLois Curfman McInnes    Input Parameters:
495c7afd0dbSLois Curfman McInnes +  snes - the SNES context
496c7afd0dbSLois Curfman McInnes -  usrP - optional user context
497c7afd0dbSLois Curfman McInnes 
49836851e7fSLois Curfman McInnes    Level: intermediate
49936851e7fSLois Curfman McInnes 
5009b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
5019b94acceSBarry Smith 
5029b94acceSBarry Smith .seealso: SNESGetApplicationContext()
5039b94acceSBarry Smith @*/
5047087cfbeSBarry Smith PetscErrorCode  SNESSetApplicationContext(SNES snes,void *usrP)
5059b94acceSBarry Smith {
5063a40ed3dSBarry Smith   PetscFunctionBegin;
5070700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5089b94acceSBarry Smith   snes->user		= usrP;
5093a40ed3dSBarry Smith   PetscFunctionReturn(0);
5109b94acceSBarry Smith }
51174679c65SBarry Smith 
5124a2ae208SSatish Balay #undef __FUNCT__
5134a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
5149b94acceSBarry Smith /*@C
5159b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
5169b94acceSBarry Smith    nonlinear solvers.
5179b94acceSBarry Smith 
518c7afd0dbSLois Curfman McInnes    Not Collective
519c7afd0dbSLois Curfman McInnes 
5209b94acceSBarry Smith    Input Parameter:
5219b94acceSBarry Smith .  snes - SNES context
5229b94acceSBarry Smith 
5239b94acceSBarry Smith    Output Parameter:
5249b94acceSBarry Smith .  usrP - user context
5259b94acceSBarry Smith 
52636851e7fSLois Curfman McInnes    Level: intermediate
52736851e7fSLois Curfman McInnes 
5289b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
5299b94acceSBarry Smith 
5309b94acceSBarry Smith .seealso: SNESSetApplicationContext()
5319b94acceSBarry Smith @*/
5327087cfbeSBarry Smith PetscErrorCode  SNESGetApplicationContext(SNES snes,void **usrP)
5339b94acceSBarry Smith {
5343a40ed3dSBarry Smith   PetscFunctionBegin;
5350700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5369b94acceSBarry Smith   *usrP = snes->user;
5373a40ed3dSBarry Smith   PetscFunctionReturn(0);
5389b94acceSBarry Smith }
53974679c65SBarry Smith 
5404a2ae208SSatish Balay #undef __FUNCT__
5414a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
5429b94acceSBarry Smith /*@
543c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
544c8228a4eSBarry Smith    at this time.
5459b94acceSBarry Smith 
546c7afd0dbSLois Curfman McInnes    Not Collective
547c7afd0dbSLois Curfman McInnes 
5489b94acceSBarry Smith    Input Parameter:
5499b94acceSBarry Smith .  snes - SNES context
5509b94acceSBarry Smith 
5519b94acceSBarry Smith    Output Parameter:
5529b94acceSBarry Smith .  iter - iteration number
5539b94acceSBarry Smith 
554c8228a4eSBarry Smith    Notes:
555c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
556c8228a4eSBarry Smith 
557c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
55808405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
55908405cd6SLois Curfman McInnes .vb
56008405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
56108405cd6SLois Curfman McInnes       if (!(it % 2)) {
56208405cd6SLois Curfman McInnes         [compute Jacobian here]
56308405cd6SLois Curfman McInnes       }
56408405cd6SLois Curfman McInnes .ve
565c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
56608405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
567c8228a4eSBarry Smith 
56836851e7fSLois Curfman McInnes    Level: intermediate
56936851e7fSLois Curfman McInnes 
5702b668275SBarry Smith .keywords: SNES, nonlinear, get, iteration, number,
5712b668275SBarry Smith 
572b850b91aSLisandro Dalcin .seealso:   SNESGetFunctionNorm(), SNESGetLinearSolveIterations()
5739b94acceSBarry Smith @*/
5747087cfbeSBarry Smith PetscErrorCode  SNESGetIterationNumber(SNES snes,PetscInt* iter)
5759b94acceSBarry Smith {
5763a40ed3dSBarry Smith   PetscFunctionBegin;
5770700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
5784482741eSBarry Smith   PetscValidIntPointer(iter,2);
5799b94acceSBarry Smith   *iter = snes->iter;
5803a40ed3dSBarry Smith   PetscFunctionReturn(0);
5819b94acceSBarry Smith }
58274679c65SBarry Smith 
5834a2ae208SSatish Balay #undef __FUNCT__
5844a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
5859b94acceSBarry Smith /*@
5869b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
5879b94acceSBarry Smith    with SNESSSetFunction().
5889b94acceSBarry Smith 
589c7afd0dbSLois Curfman McInnes    Collective on SNES
590c7afd0dbSLois Curfman McInnes 
5919b94acceSBarry Smith    Input Parameter:
5929b94acceSBarry Smith .  snes - SNES context
5939b94acceSBarry Smith 
5949b94acceSBarry Smith    Output Parameter:
5959b94acceSBarry Smith .  fnorm - 2-norm of function
5969b94acceSBarry Smith 
59736851e7fSLois Curfman McInnes    Level: intermediate
59836851e7fSLois Curfman McInnes 
5999b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
600a86d99e1SLois Curfman McInnes 
601b850b91aSLisandro Dalcin .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations()
6029b94acceSBarry Smith @*/
6037087cfbeSBarry Smith PetscErrorCode  SNESGetFunctionNorm(SNES snes,PetscReal *fnorm)
6049b94acceSBarry Smith {
6053a40ed3dSBarry Smith   PetscFunctionBegin;
6060700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6074482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
6089b94acceSBarry Smith   *fnorm = snes->norm;
6093a40ed3dSBarry Smith   PetscFunctionReturn(0);
6109b94acceSBarry Smith }
61174679c65SBarry Smith 
6124a2ae208SSatish Balay #undef __FUNCT__
613b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetNonlinearStepFailures"
6149b94acceSBarry Smith /*@
615b850b91aSLisandro Dalcin    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
6169b94acceSBarry Smith    attempted by the nonlinear solver.
6179b94acceSBarry Smith 
618c7afd0dbSLois Curfman McInnes    Not Collective
619c7afd0dbSLois Curfman McInnes 
6209b94acceSBarry Smith    Input Parameter:
6219b94acceSBarry Smith .  snes - SNES context
6229b94acceSBarry Smith 
6239b94acceSBarry Smith    Output Parameter:
6249b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
6259b94acceSBarry Smith 
626c96a6f78SLois Curfman McInnes    Notes:
627c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
628c96a6f78SLois Curfman McInnes 
62936851e7fSLois Curfman McInnes    Level: intermediate
63036851e7fSLois Curfman McInnes 
6319b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
63258ebbce7SBarry Smith 
633e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
63458ebbce7SBarry Smith           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
6359b94acceSBarry Smith @*/
6367087cfbeSBarry Smith PetscErrorCode  SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails)
6379b94acceSBarry Smith {
6383a40ed3dSBarry Smith   PetscFunctionBegin;
6390700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6404482741eSBarry Smith   PetscValidIntPointer(nfails,2);
64150ffb88aSMatthew Knepley   *nfails = snes->numFailures;
64250ffb88aSMatthew Knepley   PetscFunctionReturn(0);
64350ffb88aSMatthew Knepley }
64450ffb88aSMatthew Knepley 
64550ffb88aSMatthew Knepley #undef __FUNCT__
646b850b91aSLisandro Dalcin #define __FUNCT__ "SNESSetMaxNonlinearStepFailures"
64750ffb88aSMatthew Knepley /*@
648b850b91aSLisandro Dalcin    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
64950ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
65050ffb88aSMatthew Knepley 
65150ffb88aSMatthew Knepley    Not Collective
65250ffb88aSMatthew Knepley 
65350ffb88aSMatthew Knepley    Input Parameters:
65450ffb88aSMatthew Knepley +  snes     - SNES context
65550ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
65650ffb88aSMatthew Knepley 
65750ffb88aSMatthew Knepley    Level: intermediate
65850ffb88aSMatthew Knepley 
65950ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
66058ebbce7SBarry Smith 
661e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
66258ebbce7SBarry Smith           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
66350ffb88aSMatthew Knepley @*/
6647087cfbeSBarry Smith PetscErrorCode  SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
66550ffb88aSMatthew Knepley {
66650ffb88aSMatthew Knepley   PetscFunctionBegin;
6670700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
66850ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
66950ffb88aSMatthew Knepley   PetscFunctionReturn(0);
67050ffb88aSMatthew Knepley }
67150ffb88aSMatthew Knepley 
67250ffb88aSMatthew Knepley #undef __FUNCT__
673b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetMaxNonlinearStepFailures"
67450ffb88aSMatthew Knepley /*@
675b850b91aSLisandro Dalcin    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
67650ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
67750ffb88aSMatthew Knepley 
67850ffb88aSMatthew Knepley    Not Collective
67950ffb88aSMatthew Knepley 
68050ffb88aSMatthew Knepley    Input Parameter:
68150ffb88aSMatthew Knepley .  snes     - SNES context
68250ffb88aSMatthew Knepley 
68350ffb88aSMatthew Knepley    Output Parameter:
68450ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
68550ffb88aSMatthew Knepley 
68650ffb88aSMatthew Knepley    Level: intermediate
68750ffb88aSMatthew Knepley 
68850ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
68958ebbce7SBarry Smith 
690e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
69158ebbce7SBarry Smith           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
69258ebbce7SBarry Smith 
69350ffb88aSMatthew Knepley @*/
6947087cfbeSBarry Smith PetscErrorCode  SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
69550ffb88aSMatthew Knepley {
69650ffb88aSMatthew Knepley   PetscFunctionBegin;
6970700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6984482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
69950ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
7003a40ed3dSBarry Smith   PetscFunctionReturn(0);
7019b94acceSBarry Smith }
702a847f771SSatish Balay 
7034a2ae208SSatish Balay #undef __FUNCT__
7042541af92SBarry Smith #define __FUNCT__ "SNESGetNumberFunctionEvals"
7052541af92SBarry Smith /*@
7062541af92SBarry Smith    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
7072541af92SBarry Smith      done by SNES.
7082541af92SBarry Smith 
7092541af92SBarry Smith    Not Collective
7102541af92SBarry Smith 
7112541af92SBarry Smith    Input Parameter:
7122541af92SBarry Smith .  snes     - SNES context
7132541af92SBarry Smith 
7142541af92SBarry Smith    Output Parameter:
7152541af92SBarry Smith .  nfuncs - number of evaluations
7162541af92SBarry Smith 
7172541af92SBarry Smith    Level: intermediate
7182541af92SBarry Smith 
7192541af92SBarry Smith .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
72058ebbce7SBarry Smith 
721e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures()
7222541af92SBarry Smith @*/
7237087cfbeSBarry Smith PetscErrorCode  SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
7242541af92SBarry Smith {
7252541af92SBarry Smith   PetscFunctionBegin;
7260700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7272541af92SBarry Smith   PetscValidIntPointer(nfuncs,2);
7282541af92SBarry Smith   *nfuncs = snes->nfuncs;
7292541af92SBarry Smith   PetscFunctionReturn(0);
7302541af92SBarry Smith }
7312541af92SBarry Smith 
7322541af92SBarry Smith #undef __FUNCT__
7333d4c4710SBarry Smith #define __FUNCT__ "SNESGetLinearSolveFailures"
7343d4c4710SBarry Smith /*@
7353d4c4710SBarry Smith    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
7363d4c4710SBarry Smith    linear solvers.
7373d4c4710SBarry Smith 
7383d4c4710SBarry Smith    Not Collective
7393d4c4710SBarry Smith 
7403d4c4710SBarry Smith    Input Parameter:
7413d4c4710SBarry Smith .  snes - SNES context
7423d4c4710SBarry Smith 
7433d4c4710SBarry Smith    Output Parameter:
7443d4c4710SBarry Smith .  nfails - number of failed solves
7453d4c4710SBarry Smith 
7463d4c4710SBarry Smith    Notes:
7473d4c4710SBarry Smith    This counter is reset to zero for each successive call to SNESSolve().
7483d4c4710SBarry Smith 
7493d4c4710SBarry Smith    Level: intermediate
7503d4c4710SBarry Smith 
7513d4c4710SBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
75258ebbce7SBarry Smith 
753e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
7543d4c4710SBarry Smith @*/
7557087cfbeSBarry Smith PetscErrorCode  SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails)
7563d4c4710SBarry Smith {
7573d4c4710SBarry Smith   PetscFunctionBegin;
7580700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
7593d4c4710SBarry Smith   PetscValidIntPointer(nfails,2);
7603d4c4710SBarry Smith   *nfails = snes->numLinearSolveFailures;
7613d4c4710SBarry Smith   PetscFunctionReturn(0);
7623d4c4710SBarry Smith }
7633d4c4710SBarry Smith 
7643d4c4710SBarry Smith #undef __FUNCT__
7653d4c4710SBarry Smith #define __FUNCT__ "SNESSetMaxLinearSolveFailures"
7663d4c4710SBarry Smith /*@
7673d4c4710SBarry Smith    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
7683d4c4710SBarry Smith    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
7693d4c4710SBarry Smith 
7703f9fe445SBarry Smith    Logically Collective on SNES
7713d4c4710SBarry Smith 
7723d4c4710SBarry Smith    Input Parameters:
7733d4c4710SBarry Smith +  snes     - SNES context
7743d4c4710SBarry Smith -  maxFails - maximum allowed linear solve failures
7753d4c4710SBarry Smith 
7763d4c4710SBarry Smith    Level: intermediate
7773d4c4710SBarry Smith 
778a6796414SBarry Smith    Notes: By default this is 0; that is SNES returns on the first failed linear solve
7793d4c4710SBarry Smith 
7803d4c4710SBarry Smith .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
7813d4c4710SBarry Smith 
78258ebbce7SBarry Smith .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
7833d4c4710SBarry Smith @*/
7847087cfbeSBarry Smith PetscErrorCode  SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
7853d4c4710SBarry Smith {
7863d4c4710SBarry Smith   PetscFunctionBegin;
7870700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
788c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(snes,maxFails,2);
7893d4c4710SBarry Smith   snes->maxLinearSolveFailures = maxFails;
7903d4c4710SBarry Smith   PetscFunctionReturn(0);
7913d4c4710SBarry Smith }
7923d4c4710SBarry Smith 
7933d4c4710SBarry Smith #undef __FUNCT__
7943d4c4710SBarry Smith #define __FUNCT__ "SNESGetMaxLinearSolveFailures"
7953d4c4710SBarry Smith /*@
7963d4c4710SBarry Smith    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
7973d4c4710SBarry Smith      are allowed before SNES terminates
7983d4c4710SBarry Smith 
7993d4c4710SBarry Smith    Not Collective
8003d4c4710SBarry Smith 
8013d4c4710SBarry Smith    Input Parameter:
8023d4c4710SBarry Smith .  snes     - SNES context
8033d4c4710SBarry Smith 
8043d4c4710SBarry Smith    Output Parameter:
8053d4c4710SBarry Smith .  maxFails - maximum of unsuccessful solves allowed
8063d4c4710SBarry Smith 
8073d4c4710SBarry Smith    Level: intermediate
8083d4c4710SBarry Smith 
8093d4c4710SBarry Smith    Notes: By default this is 1; that is SNES returns on the first failed linear solve
8103d4c4710SBarry Smith 
8113d4c4710SBarry Smith .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
8123d4c4710SBarry Smith 
813e1c61ce8SBarry Smith .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
8143d4c4710SBarry Smith @*/
8157087cfbeSBarry Smith PetscErrorCode  SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
8163d4c4710SBarry Smith {
8173d4c4710SBarry Smith   PetscFunctionBegin;
8180700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
8193d4c4710SBarry Smith   PetscValidIntPointer(maxFails,2);
8203d4c4710SBarry Smith   *maxFails = snes->maxLinearSolveFailures;
8213d4c4710SBarry Smith   PetscFunctionReturn(0);
8223d4c4710SBarry Smith }
8233d4c4710SBarry Smith 
8243d4c4710SBarry Smith #undef __FUNCT__
825b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetLinearSolveIterations"
826c96a6f78SLois Curfman McInnes /*@
827b850b91aSLisandro Dalcin    SNESGetLinearSolveIterations - Gets the total number of linear iterations
828c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
829c96a6f78SLois Curfman McInnes 
830c7afd0dbSLois Curfman McInnes    Not Collective
831c7afd0dbSLois Curfman McInnes 
832c96a6f78SLois Curfman McInnes    Input Parameter:
833c96a6f78SLois Curfman McInnes .  snes - SNES context
834c96a6f78SLois Curfman McInnes 
835c96a6f78SLois Curfman McInnes    Output Parameter:
836c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
837c96a6f78SLois Curfman McInnes 
838c96a6f78SLois Curfman McInnes    Notes:
839c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
840c96a6f78SLois Curfman McInnes 
84136851e7fSLois Curfman McInnes    Level: intermediate
84236851e7fSLois Curfman McInnes 
843c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
8442b668275SBarry Smith 
8458c7482ecSBarry Smith .seealso:  SNESGetIterationNumber(), SNESGetFunctionNorm(), SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures()
846c96a6f78SLois Curfman McInnes @*/
8477087cfbeSBarry Smith PetscErrorCode  SNESGetLinearSolveIterations(SNES snes,PetscInt* lits)
848c96a6f78SLois Curfman McInnes {
8493a40ed3dSBarry Smith   PetscFunctionBegin;
8500700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
8514482741eSBarry Smith   PetscValidIntPointer(lits,2);
852c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
8533a40ed3dSBarry Smith   PetscFunctionReturn(0);
854c96a6f78SLois Curfman McInnes }
855c96a6f78SLois Curfman McInnes 
8564a2ae208SSatish Balay #undef __FUNCT__
85794b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
85852baeb72SSatish Balay /*@
85994b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
8609b94acceSBarry Smith 
86194b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
862c7afd0dbSLois Curfman McInnes 
8639b94acceSBarry Smith    Input Parameter:
8649b94acceSBarry Smith .  snes - the SNES context
8659b94acceSBarry Smith 
8669b94acceSBarry Smith    Output Parameter:
86794b7f48cSBarry Smith .  ksp - the KSP context
8689b94acceSBarry Smith 
8699b94acceSBarry Smith    Notes:
87094b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
8719b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
8722999313aSBarry Smith    PC contexts as well.
8739b94acceSBarry Smith 
87436851e7fSLois Curfman McInnes    Level: beginner
87536851e7fSLois Curfman McInnes 
87694b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
8779b94acceSBarry Smith 
8782999313aSBarry Smith .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
8799b94acceSBarry Smith @*/
8807087cfbeSBarry Smith PetscErrorCode  SNESGetKSP(SNES snes,KSP *ksp)
8819b94acceSBarry Smith {
8821cee3971SBarry Smith   PetscErrorCode ierr;
8831cee3971SBarry Smith 
8843a40ed3dSBarry Smith   PetscFunctionBegin;
8850700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
8864482741eSBarry Smith   PetscValidPointer(ksp,2);
8871cee3971SBarry Smith 
8881cee3971SBarry Smith   if (!snes->ksp) {
8891cee3971SBarry Smith     ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr);
8901cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr);
8911cee3971SBarry Smith     ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
8921cee3971SBarry Smith   }
89394b7f48cSBarry Smith   *ksp = snes->ksp;
8943a40ed3dSBarry Smith   PetscFunctionReturn(0);
8959b94acceSBarry Smith }
89682bf6240SBarry Smith 
8974a2ae208SSatish Balay #undef __FUNCT__
8982999313aSBarry Smith #define __FUNCT__ "SNESSetKSP"
8992999313aSBarry Smith /*@
9002999313aSBarry Smith    SNESSetKSP - Sets a KSP context for the SNES object to use
9012999313aSBarry Smith 
9022999313aSBarry Smith    Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
9032999313aSBarry Smith 
9042999313aSBarry Smith    Input Parameters:
9052999313aSBarry Smith +  snes - the SNES context
9062999313aSBarry Smith -  ksp - the KSP context
9072999313aSBarry Smith 
9082999313aSBarry Smith    Notes:
9092999313aSBarry Smith    The SNES object already has its KSP object, you can obtain with SNESGetKSP()
9102999313aSBarry Smith    so this routine is rarely needed.
9112999313aSBarry Smith 
9122999313aSBarry Smith    The KSP object that is already in the SNES object has its reference count
9132999313aSBarry Smith    decreased by one.
9142999313aSBarry Smith 
9152999313aSBarry Smith    Level: developer
9162999313aSBarry Smith 
9172999313aSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
9182999313aSBarry Smith 
9192999313aSBarry Smith .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
9202999313aSBarry Smith @*/
9217087cfbeSBarry Smith PetscErrorCode  SNESSetKSP(SNES snes,KSP ksp)
9222999313aSBarry Smith {
9232999313aSBarry Smith   PetscErrorCode ierr;
9242999313aSBarry Smith 
9252999313aSBarry Smith   PetscFunctionBegin;
9260700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
9270700a824SBarry Smith   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
9282999313aSBarry Smith   PetscCheckSameComm(snes,1,ksp,2);
9297dcf0eaaSdalcinl   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
930906ed7ccSBarry Smith   if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);}
9312999313aSBarry Smith   snes->ksp = ksp;
9322999313aSBarry Smith   PetscFunctionReturn(0);
9332999313aSBarry Smith }
9342999313aSBarry Smith 
9357adad957SLisandro Dalcin #if 0
9362999313aSBarry Smith #undef __FUNCT__
9374a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
9386849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
939e24b481bSBarry Smith {
940e24b481bSBarry Smith   PetscFunctionBegin;
941e24b481bSBarry Smith   PetscFunctionReturn(0);
942e24b481bSBarry Smith }
9437adad957SLisandro Dalcin #endif
944e24b481bSBarry Smith 
9459b94acceSBarry Smith /* -----------------------------------------------------------*/
9464a2ae208SSatish Balay #undef __FUNCT__
9474a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
94852baeb72SSatish Balay /*@
9499b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
9509b94acceSBarry Smith 
951c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
952c7afd0dbSLois Curfman McInnes 
953c7afd0dbSLois Curfman McInnes    Input Parameters:
954906ed7ccSBarry Smith .  comm - MPI communicator
9559b94acceSBarry Smith 
9569b94acceSBarry Smith    Output Parameter:
9579b94acceSBarry Smith .  outsnes - the new SNES context
9589b94acceSBarry Smith 
959c7afd0dbSLois Curfman McInnes    Options Database Keys:
960c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
961c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
962c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
963c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
964c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
965c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
966c1f60f51SBarry Smith 
96736851e7fSLois Curfman McInnes    Level: beginner
96836851e7fSLois Curfman McInnes 
9699b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
9709b94acceSBarry Smith 
971a8054027SBarry Smith .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
972a8054027SBarry Smith 
9739b94acceSBarry Smith @*/
9747087cfbeSBarry Smith PetscErrorCode  SNESCreate(MPI_Comm comm,SNES *outsnes)
9759b94acceSBarry Smith {
976dfbe8321SBarry Smith   PetscErrorCode      ierr;
9779b94acceSBarry Smith   SNES                snes;
978fa9f3622SBarry Smith   SNESKSPEW           *kctx;
97937fcc0dbSBarry Smith 
9803a40ed3dSBarry Smith   PetscFunctionBegin;
981ed1caa07SMatthew Knepley   PetscValidPointer(outsnes,2);
9828ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
9838ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
9848ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
9858ba1e511SMatthew Knepley #endif
9868ba1e511SMatthew Knepley 
9870700a824SBarry Smith   ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_CLASSID,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
9887adad957SLisandro Dalcin 
98985385478SLisandro Dalcin   snes->ops->converged    = SNESDefaultConverged;
9909b94acceSBarry Smith   snes->max_its           = 50;
9919750a799SBarry Smith   snes->max_funcs	  = 10000;
9929b94acceSBarry Smith   snes->norm		  = 0.0;
993b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
994b4874afaSBarry Smith   snes->ttol              = 0.0;
99570441072SBarry Smith   snes->abstol		  = 1.e-50;
9969b94acceSBarry Smith   snes->xtol		  = 1.e-8;
9974b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
9989b94acceSBarry Smith   snes->nfuncs            = 0;
99950ffb88aSMatthew Knepley   snes->numFailures       = 0;
100050ffb88aSMatthew Knepley   snes->maxFailures       = 1;
10017a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
1002e35cf81dSBarry Smith   snes->lagjacobian       = 1;
1003a8054027SBarry Smith   snes->lagpreconditioner = 1;
1004639f9d9dSBarry Smith   snes->numbermonitors    = 0;
10059b94acceSBarry Smith   snes->data              = 0;
10064dc4c822SBarry Smith   snes->setupcalled       = PETSC_FALSE;
1007186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
10086f24a144SLois Curfman McInnes   snes->vwork             = 0;
10096f24a144SLois Curfman McInnes   snes->nwork             = 0;
1010758f92a0SBarry Smith   snes->conv_hist_len     = 0;
1011758f92a0SBarry Smith   snes->conv_hist_max     = 0;
1012758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
1013758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
1014758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
1015184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
10169b94acceSBarry Smith 
10173d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
10183d4c4710SBarry Smith   snes->maxLinearSolveFailures = 1;
10193d4c4710SBarry Smith 
10209b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
102138f2d2fdSLisandro Dalcin   ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr);
10229b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
10239b94acceSBarry Smith   kctx->version     = 2;
10249b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
10259b94acceSBarry Smith                              this was too large for some test cases */
102675567043SBarry Smith   kctx->rtol_last   = 0.0;
10279b94acceSBarry Smith   kctx->rtol_max    = .9;
10289b94acceSBarry Smith   kctx->gamma       = 1.0;
102971f87433Sdalcinl   kctx->alpha       = .5*(1.0 + sqrt(5.0));
103071f87433Sdalcinl   kctx->alpha2      = kctx->alpha;
10319b94acceSBarry Smith   kctx->threshold   = .1;
103275567043SBarry Smith   kctx->lresid_last = 0.0;
103375567043SBarry Smith   kctx->norm_last   = 0.0;
10349b94acceSBarry Smith 
10359b94acceSBarry Smith   *outsnes = snes;
10363a40ed3dSBarry Smith   PetscFunctionReturn(0);
10379b94acceSBarry Smith }
10389b94acceSBarry Smith 
10394a2ae208SSatish Balay #undef __FUNCT__
10404a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
10419b94acceSBarry Smith /*@C
10429b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
10439b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
10449b94acceSBarry Smith    equations.
10459b94acceSBarry Smith 
10463f9fe445SBarry Smith    Logically Collective on SNES
1047fee21e36SBarry Smith 
1048c7afd0dbSLois Curfman McInnes    Input Parameters:
1049c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1050c7afd0dbSLois Curfman McInnes .  r - vector to store function value
1051de044059SHong Zhang .  func - function evaluation routine
1052c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
1053c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
10549b94acceSBarry Smith 
1055c7afd0dbSLois Curfman McInnes    Calling sequence of func:
10568d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
1057c7afd0dbSLois Curfman McInnes 
1058313e4042SLois Curfman McInnes .  f - function vector
1059c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
10609b94acceSBarry Smith 
10619b94acceSBarry Smith    Notes:
10629b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
10639b94acceSBarry Smith $      f'(x) x = -f(x),
1064c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
10659b94acceSBarry Smith 
106636851e7fSLois Curfman McInnes    Level: beginner
106736851e7fSLois Curfman McInnes 
10689b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
10699b94acceSBarry Smith 
1070a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
10719b94acceSBarry Smith @*/
10727087cfbeSBarry Smith PetscErrorCode  SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
10739b94acceSBarry Smith {
107485385478SLisandro Dalcin   PetscErrorCode ierr;
10753a40ed3dSBarry Smith   PetscFunctionBegin;
10760700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1077ef8dffc7SBarry Smith   if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2);
1078ef8dffc7SBarry Smith   if (r) PetscCheckSameComm(snes,1,r,2);
1079ef8dffc7SBarry Smith   if (!r && snes->dm) {
1080ef8dffc7SBarry Smith     ierr = DMCreateGlobalVector(snes->dm,&r);CHKERRQ(ierr);
1081ef8dffc7SBarry Smith   } else {
108285385478SLisandro Dalcin     ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr);
1083ef8dffc7SBarry Smith   }
108485385478SLisandro Dalcin   if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); }
1085e7788613SBarry Smith   snes->ops->computefunction = func;
108685385478SLisandro Dalcin   snes->vec_func             = r;
10879b94acceSBarry Smith   snes->funP                 = ctx;
10883a40ed3dSBarry Smith   PetscFunctionReturn(0);
10899b94acceSBarry Smith }
10909b94acceSBarry Smith 
10913ab0aad5SBarry Smith /* --------------------------------------------------------------- */
10923ab0aad5SBarry Smith #undef __FUNCT__
10931096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs"
10941096aae1SMatthew Knepley /*@C
10951096aae1SMatthew Knepley    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
10961096aae1SMatthew Knepley    it assumes a zero right hand side.
10971096aae1SMatthew Knepley 
10983f9fe445SBarry Smith    Logically Collective on SNES
10991096aae1SMatthew Knepley 
11001096aae1SMatthew Knepley    Input Parameter:
11011096aae1SMatthew Knepley .  snes - the SNES context
11021096aae1SMatthew Knepley 
11031096aae1SMatthew Knepley    Output Parameter:
1104bc08b0f1SBarry Smith .  rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null
11051096aae1SMatthew Knepley 
11061096aae1SMatthew Knepley    Level: intermediate
11071096aae1SMatthew Knepley 
11081096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side
11091096aae1SMatthew Knepley 
111085385478SLisandro Dalcin .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
11111096aae1SMatthew Knepley @*/
11127087cfbeSBarry Smith PetscErrorCode  SNESGetRhs(SNES snes,Vec *rhs)
11131096aae1SMatthew Knepley {
11141096aae1SMatthew Knepley   PetscFunctionBegin;
11150700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
11161096aae1SMatthew Knepley   PetscValidPointer(rhs,2);
111785385478SLisandro Dalcin   *rhs = snes->vec_rhs;
11181096aae1SMatthew Knepley   PetscFunctionReturn(0);
11191096aae1SMatthew Knepley }
11201096aae1SMatthew Knepley 
11211096aae1SMatthew Knepley #undef __FUNCT__
11224a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
11239b94acceSBarry Smith /*@
112436851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
11259b94acceSBarry Smith                          SNESSetFunction().
11269b94acceSBarry Smith 
1127c7afd0dbSLois Curfman McInnes    Collective on SNES
1128c7afd0dbSLois Curfman McInnes 
11299b94acceSBarry Smith    Input Parameters:
1130c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1131c7afd0dbSLois Curfman McInnes -  x - input vector
11329b94acceSBarry Smith 
11339b94acceSBarry Smith    Output Parameter:
11343638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
11359b94acceSBarry Smith 
11361bffabb2SLois Curfman McInnes    Notes:
113736851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
113836851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
113936851e7fSLois Curfman McInnes    themselves.
114036851e7fSLois Curfman McInnes 
114136851e7fSLois Curfman McInnes    Level: developer
114236851e7fSLois Curfman McInnes 
11439b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
11449b94acceSBarry Smith 
1145a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
11469b94acceSBarry Smith @*/
11477087cfbeSBarry Smith PetscErrorCode  SNESComputeFunction(SNES snes,Vec x,Vec y)
11489b94acceSBarry Smith {
1149dfbe8321SBarry Smith   PetscErrorCode ierr;
11509b94acceSBarry Smith 
11513a40ed3dSBarry Smith   PetscFunctionBegin;
11520700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
11530700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
11540700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
1155c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
1156c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
1157184914b5SBarry Smith 
1158d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
1159e7788613SBarry Smith   if (snes->ops->computefunction) {
1160d64ed03dSBarry Smith     PetscStackPush("SNES user function");
1161*39d508bbSBarry Smith     ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
1162d64ed03dSBarry Smith     PetscStackPop;
116385385478SLisandro Dalcin   } else if (snes->vec_rhs) {
11641096aae1SMatthew Knepley     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
116517186662SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
116685385478SLisandro Dalcin   if (snes->vec_rhs) {
116785385478SLisandro Dalcin     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
11683ab0aad5SBarry Smith   }
1169ae3c334cSLois Curfman McInnes   snes->nfuncs++;
1170d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
11713a40ed3dSBarry Smith   PetscFunctionReturn(0);
11729b94acceSBarry Smith }
11739b94acceSBarry Smith 
11744a2ae208SSatish Balay #undef __FUNCT__
11754a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
117662fef451SLois Curfman McInnes /*@
117762fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
117862fef451SLois Curfman McInnes    set with SNESSetJacobian().
117962fef451SLois Curfman McInnes 
1180c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1181c7afd0dbSLois Curfman McInnes 
118262fef451SLois Curfman McInnes    Input Parameters:
1183c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1184c7afd0dbSLois Curfman McInnes -  x - input vector
118562fef451SLois Curfman McInnes 
118662fef451SLois Curfman McInnes    Output Parameters:
1187c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
118862fef451SLois Curfman McInnes .  B - optional preconditioning matrix
11892b668275SBarry Smith -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1190fee21e36SBarry Smith 
1191e35cf81dSBarry Smith   Options Database Keys:
1192e35cf81dSBarry Smith +    -snes_lag_preconditioner <lag>
1193e35cf81dSBarry Smith -    -snes_lag_jacobian <lag>
1194e35cf81dSBarry Smith 
119562fef451SLois Curfman McInnes    Notes:
119662fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
119762fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
119862fef451SLois Curfman McInnes 
119994b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
1200dc5a77f8SLois Curfman McInnes    flag parameter.
120162fef451SLois Curfman McInnes 
120236851e7fSLois Curfman McInnes    Level: developer
120336851e7fSLois Curfman McInnes 
120462fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
120562fef451SLois Curfman McInnes 
1206e35cf81dSBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian()
120762fef451SLois Curfman McInnes @*/
12087087cfbeSBarry Smith PetscErrorCode  SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
12099b94acceSBarry Smith {
1210dfbe8321SBarry Smith   PetscErrorCode ierr;
1211ace3abfcSBarry Smith   PetscBool      flag;
12123a40ed3dSBarry Smith 
12133a40ed3dSBarry Smith   PetscFunctionBegin;
12140700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
12150700a824SBarry Smith   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
12164482741eSBarry Smith   PetscValidPointer(flg,5);
1217c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
1218e7788613SBarry Smith   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1219ebd3b9afSBarry Smith 
1220ebd3b9afSBarry Smith   /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */
1221ebd3b9afSBarry Smith 
1222fe3ffe1eSBarry Smith   if (snes->lagjacobian == -2) {
1223fe3ffe1eSBarry Smith     snes->lagjacobian = -1;
1224fe3ffe1eSBarry Smith     ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr);
1225fe3ffe1eSBarry Smith   } else if (snes->lagjacobian == -1) {
1226e35cf81dSBarry Smith     *flg = SAME_PRECONDITIONER;
1227e35cf81dSBarry Smith     ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr);
1228ebd3b9afSBarry Smith     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1229ebd3b9afSBarry Smith     if (flag) {
1230ebd3b9afSBarry Smith       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1231ebd3b9afSBarry Smith       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1232ebd3b9afSBarry Smith     }
1233e35cf81dSBarry Smith     PetscFunctionReturn(0);
1234e35cf81dSBarry Smith   } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) {
1235e35cf81dSBarry Smith     *flg = SAME_PRECONDITIONER;
1236e35cf81dSBarry Smith     ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr);
1237ebd3b9afSBarry Smith     ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr);
1238ebd3b9afSBarry Smith     if (flag) {
1239ebd3b9afSBarry Smith       ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1240ebd3b9afSBarry Smith       ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1241ebd3b9afSBarry Smith     }
1242e35cf81dSBarry Smith     PetscFunctionReturn(0);
1243e35cf81dSBarry Smith   }
1244e35cf81dSBarry Smith 
1245c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
1246e35cf81dSBarry Smith   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1247d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
1248e7788613SBarry Smith   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1249d64ed03dSBarry Smith   PetscStackPop;
1250d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1251a8054027SBarry Smith 
12523b4f5425SBarry Smith   if (snes->lagpreconditioner == -2) {
12533b4f5425SBarry Smith     ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr);
12543b4f5425SBarry Smith     snes->lagpreconditioner = -1;
12553b4f5425SBarry Smith   } else if (snes->lagpreconditioner == -1) {
1256a8054027SBarry Smith     *flg = SAME_PRECONDITIONER;
1257a8054027SBarry Smith     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1258a8054027SBarry Smith   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1259a8054027SBarry Smith     *flg = SAME_PRECONDITIONER;
1260a8054027SBarry Smith     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1261a8054027SBarry Smith   }
1262a8054027SBarry Smith 
12636d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
12640700a824SBarry Smith   /* PetscValidHeaderSpecific(*A,MAT_CLASSID,3);
12650700a824SBarry Smith     PetscValidHeaderSpecific(*B,MAT_CLASSID,4);   */
12663a40ed3dSBarry Smith   PetscFunctionReturn(0);
12679b94acceSBarry Smith }
12689b94acceSBarry Smith 
12694a2ae208SSatish Balay #undef __FUNCT__
12704a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
12719b94acceSBarry Smith /*@C
12729b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1273044dda88SLois Curfman McInnes    location to store the matrix.
12749b94acceSBarry Smith 
12753f9fe445SBarry Smith    Logically Collective on SNES and Mat
1276c7afd0dbSLois Curfman McInnes 
12779b94acceSBarry Smith    Input Parameters:
1278c7afd0dbSLois Curfman McInnes +  snes - the SNES context
12799b94acceSBarry Smith .  A - Jacobian matrix
12809b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
12819b94acceSBarry Smith .  func - Jacobian evaluation routine
1282c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
12832cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
12849b94acceSBarry Smith 
12859b94acceSBarry Smith    Calling sequence of func:
12868d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
12879b94acceSBarry Smith 
1288c7afd0dbSLois Curfman McInnes +  x - input vector
12899b94acceSBarry Smith .  A - Jacobian matrix
12909b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1291ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
12922b668275SBarry Smith    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1293c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
12949b94acceSBarry Smith 
12959b94acceSBarry Smith    Notes:
129694b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
12972cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1298ac21db08SLois Curfman McInnes 
1299ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
13009b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
13019b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
13029b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
13039b94acceSBarry Smith    throughout the global iterations.
13049b94acceSBarry Smith 
130516913363SBarry Smith    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
130616913363SBarry Smith    each matrix.
130716913363SBarry Smith 
130836851e7fSLois Curfman McInnes    Level: beginner
130936851e7fSLois Curfman McInnes 
13109b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
13119b94acceSBarry Smith 
13123ec795f1SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
13139b94acceSBarry Smith @*/
13147087cfbeSBarry Smith PetscErrorCode  SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
13159b94acceSBarry Smith {
1316dfbe8321SBarry Smith   PetscErrorCode ierr;
13173a7fca6bSBarry Smith 
13183a40ed3dSBarry Smith   PetscFunctionBegin;
13190700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
13200700a824SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2);
13210700a824SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3);
1322c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
132306975374SJed Brown   if (B) PetscCheckSameComm(snes,1,B,3);
1324e7788613SBarry Smith   if (func) snes->ops->computejacobian = func;
13253a7fca6bSBarry Smith   if (ctx)  snes->jacP                 = ctx;
13263a7fca6bSBarry Smith   if (A) {
13277dcf0eaaSdalcinl     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
13283a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
13299b94acceSBarry Smith     snes->jacobian = A;
13303a7fca6bSBarry Smith   }
13313a7fca6bSBarry Smith   if (B) {
13327dcf0eaaSdalcinl     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
13333a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
13349b94acceSBarry Smith     snes->jacobian_pre = B;
13353a7fca6bSBarry Smith   }
13363a40ed3dSBarry Smith   PetscFunctionReturn(0);
13379b94acceSBarry Smith }
133862fef451SLois Curfman McInnes 
13394a2ae208SSatish Balay #undef __FUNCT__
13404a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
1341c2aafc4cSSatish Balay /*@C
1342b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1343b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1344b4fd4287SBarry Smith 
1345c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1346c7afd0dbSLois Curfman McInnes 
1347b4fd4287SBarry Smith    Input Parameter:
1348b4fd4287SBarry Smith .  snes - the nonlinear solver context
1349b4fd4287SBarry Smith 
1350b4fd4287SBarry Smith    Output Parameters:
1351c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1352b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
135370e92668SMatthew Knepley .  func - location to put Jacobian function (or PETSC_NULL)
135470e92668SMatthew Knepley -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1355fee21e36SBarry Smith 
135636851e7fSLois Curfman McInnes    Level: advanced
135736851e7fSLois Curfman McInnes 
1358b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1359b4fd4287SBarry Smith @*/
13607087cfbeSBarry Smith PetscErrorCode  SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1361b4fd4287SBarry Smith {
13623a40ed3dSBarry Smith   PetscFunctionBegin;
13630700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1364b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
1365b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
1366e7788613SBarry Smith   if (func) *func = snes->ops->computejacobian;
136770e92668SMatthew Knepley   if (ctx)  *ctx  = snes->jacP;
13683a40ed3dSBarry Smith   PetscFunctionReturn(0);
1369b4fd4287SBarry Smith }
1370b4fd4287SBarry Smith 
13719b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
13729b94acceSBarry Smith 
13734a2ae208SSatish Balay #undef __FUNCT__
13744a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
13759b94acceSBarry Smith /*@
13769b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1377272ac6f2SLois Curfman McInnes    of a nonlinear solver.
13789b94acceSBarry Smith 
1379fee21e36SBarry Smith    Collective on SNES
1380fee21e36SBarry Smith 
1381c7afd0dbSLois Curfman McInnes    Input Parameters:
138270e92668SMatthew Knepley .  snes - the SNES context
1383c7afd0dbSLois Curfman McInnes 
1384272ac6f2SLois Curfman McInnes    Notes:
1385272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1386272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1387272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1388272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1389272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1390272ac6f2SLois Curfman McInnes 
139136851e7fSLois Curfman McInnes    Level: advanced
139236851e7fSLois Curfman McInnes 
13939b94acceSBarry Smith .keywords: SNES, nonlinear, setup
13949b94acceSBarry Smith 
13959b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
13969b94acceSBarry Smith @*/
13977087cfbeSBarry Smith PetscErrorCode  SNESSetUp(SNES snes)
13989b94acceSBarry Smith {
1399dfbe8321SBarry Smith   PetscErrorCode ierr;
14003a40ed3dSBarry Smith 
14013a40ed3dSBarry Smith   PetscFunctionBegin;
14020700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
14034dc4c822SBarry Smith   if (snes->setupcalled) PetscFunctionReturn(0);
14049b94acceSBarry Smith 
14057adad957SLisandro Dalcin   if (!((PetscObject)snes)->type_name) {
140685385478SLisandro Dalcin     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
140785385478SLisandro Dalcin   }
140885385478SLisandro Dalcin 
140917186662SBarry Smith   if (!snes->vec_func && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
141017186662SBarry Smith   if (!snes->ops->computefunction && !snes->vec_rhs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
141117186662SBarry Smith   if (snes->vec_func == snes->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1412ef8dffc7SBarry Smith   if (!snes->ops->computejacobian && snes->dm) {
1413ef8dffc7SBarry Smith     Mat           J;
1414ef8dffc7SBarry Smith     ISColoring    coloring;
1415ef8dffc7SBarry Smith     MatFDColoring fd;
1416a703fe33SLois Curfman McInnes 
1417ef8dffc7SBarry Smith     ierr = DMGetMatrix(snes->dm,MATAIJ,&J);CHKERRQ(ierr);
1418ef8dffc7SBarry Smith     ierr = DMGetColoring(snes->dm,IS_COLORING_GHOSTED,MATAIJ,&coloring);CHKERRQ(ierr);
1419ef8dffc7SBarry Smith     ierr = MatFDColoringCreate(J,coloring,&fd);CHKERRQ(ierr);
1420ef8dffc7SBarry Smith     ierr = MatFDColoringSetFunction(fd,(PetscErrorCode (*)(void))snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
1421ef8dffc7SBarry Smith     ierr = SNESSetJacobian(snes,J,J,SNESDefaultComputeJacobianColor,fd);CHKERRQ(ierr);
1422ef8dffc7SBarry Smith     ierr = ISColoringDestroy(coloring);CHKERRQ(ierr);
1423ef8dffc7SBarry Smith   }
1424b710008aSBarry Smith   if (!snes->ksp) {ierr = SNESGetKSP(snes, &snes->ksp);CHKERRQ(ierr);}
1425b710008aSBarry Smith 
1426410397dcSLisandro Dalcin   if (snes->ops->setup) {
1427410397dcSLisandro Dalcin     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1428410397dcSLisandro Dalcin   }
14297aaed0d8SBarry Smith   snes->setupcalled = PETSC_TRUE;
14303a40ed3dSBarry Smith   PetscFunctionReturn(0);
14319b94acceSBarry Smith }
14329b94acceSBarry Smith 
14334a2ae208SSatish Balay #undef __FUNCT__
14344a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
143552baeb72SSatish Balay /*@
14369b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
14379b94acceSBarry Smith    with SNESCreate().
14389b94acceSBarry Smith 
1439c7afd0dbSLois Curfman McInnes    Collective on SNES
1440c7afd0dbSLois Curfman McInnes 
14419b94acceSBarry Smith    Input Parameter:
14429b94acceSBarry Smith .  snes - the SNES context
14439b94acceSBarry Smith 
144436851e7fSLois Curfman McInnes    Level: beginner
144536851e7fSLois Curfman McInnes 
14469b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
14479b94acceSBarry Smith 
144863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
14499b94acceSBarry Smith @*/
14507087cfbeSBarry Smith PetscErrorCode  SNESDestroy(SNES snes)
14519b94acceSBarry Smith {
14526849ba73SBarry Smith   PetscErrorCode ierr;
14533a40ed3dSBarry Smith 
14543a40ed3dSBarry Smith   PetscFunctionBegin;
14550700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
14567adad957SLisandro Dalcin   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1457d4bb536fSBarry Smith 
1458be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
14590f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
14606c699258SBarry Smith 
14616c699258SBarry Smith   if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);}
1462e7788613SBarry Smith   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
146385385478SLisandro Dalcin 
146485385478SLisandro Dalcin   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
146585385478SLisandro Dalcin   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
146685385478SLisandro Dalcin   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
14673a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
14683a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
14691cee3971SBarry Smith   if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);}
147085385478SLisandro Dalcin   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1471522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1472a6570f20SBarry Smith   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
14737f7931b9SBarry Smith   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1474a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
14753a40ed3dSBarry Smith   PetscFunctionReturn(0);
14769b94acceSBarry Smith }
14779b94acceSBarry Smith 
14789b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
14799b94acceSBarry Smith 
14804a2ae208SSatish Balay #undef __FUNCT__
1481a8054027SBarry Smith #define __FUNCT__ "SNESSetLagPreconditioner"
1482a8054027SBarry Smith /*@
1483a8054027SBarry Smith    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1484a8054027SBarry Smith 
14853f9fe445SBarry Smith    Logically Collective on SNES
1486a8054027SBarry Smith 
1487a8054027SBarry Smith    Input Parameters:
1488a8054027SBarry Smith +  snes - the SNES context
1489a8054027SBarry Smith -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
14903b4f5425SBarry Smith          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1491a8054027SBarry Smith 
1492a8054027SBarry Smith    Options Database Keys:
1493a8054027SBarry Smith .    -snes_lag_preconditioner <lag>
1494a8054027SBarry Smith 
1495a8054027SBarry Smith    Notes:
1496a8054027SBarry Smith    The default is 1
1497a8054027SBarry Smith    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1498a8054027SBarry Smith    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1499a8054027SBarry Smith 
1500a8054027SBarry Smith    Level: intermediate
1501a8054027SBarry Smith 
1502a8054027SBarry Smith .keywords: SNES, nonlinear, set, convergence, tolerances
1503a8054027SBarry Smith 
1504e35cf81dSBarry Smith .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian()
1505a8054027SBarry Smith 
1506a8054027SBarry Smith @*/
15077087cfbeSBarry Smith PetscErrorCode  SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1508a8054027SBarry Smith {
1509a8054027SBarry Smith   PetscFunctionBegin;
15100700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1511e32f2f54SBarry Smith   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1512e32f2f54SBarry Smith   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1513c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(snes,lag,2);
1514a8054027SBarry Smith   snes->lagpreconditioner = lag;
1515a8054027SBarry Smith   PetscFunctionReturn(0);
1516a8054027SBarry Smith }
1517a8054027SBarry Smith 
1518a8054027SBarry Smith #undef __FUNCT__
1519a8054027SBarry Smith #define __FUNCT__ "SNESGetLagPreconditioner"
1520a8054027SBarry Smith /*@
1521a8054027SBarry Smith    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1522a8054027SBarry Smith 
15233f9fe445SBarry Smith    Not Collective
1524a8054027SBarry Smith 
1525a8054027SBarry Smith    Input Parameter:
1526a8054027SBarry Smith .  snes - the SNES context
1527a8054027SBarry Smith 
1528a8054027SBarry Smith    Output Parameter:
1529a8054027SBarry Smith .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
15303b4f5425SBarry Smith          the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that
1531a8054027SBarry Smith 
1532a8054027SBarry Smith    Options Database Keys:
1533a8054027SBarry Smith .    -snes_lag_preconditioner <lag>
1534a8054027SBarry Smith 
1535a8054027SBarry Smith    Notes:
1536a8054027SBarry Smith    The default is 1
1537a8054027SBarry Smith    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1538a8054027SBarry Smith 
1539a8054027SBarry Smith    Level: intermediate
1540a8054027SBarry Smith 
1541a8054027SBarry Smith .keywords: SNES, nonlinear, set, convergence, tolerances
1542a8054027SBarry Smith 
1543a8054027SBarry Smith .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1544a8054027SBarry Smith 
1545a8054027SBarry Smith @*/
15467087cfbeSBarry Smith PetscErrorCode  SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1547a8054027SBarry Smith {
1548a8054027SBarry Smith   PetscFunctionBegin;
15490700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1550a8054027SBarry Smith   *lag = snes->lagpreconditioner;
1551a8054027SBarry Smith   PetscFunctionReturn(0);
1552a8054027SBarry Smith }
1553a8054027SBarry Smith 
1554a8054027SBarry Smith #undef __FUNCT__
1555e35cf81dSBarry Smith #define __FUNCT__ "SNESSetLagJacobian"
1556e35cf81dSBarry Smith /*@
1557e35cf81dSBarry Smith    SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how
1558e35cf81dSBarry Smith      often the preconditioner is rebuilt.
1559e35cf81dSBarry Smith 
15603f9fe445SBarry Smith    Logically Collective on SNES
1561e35cf81dSBarry Smith 
1562e35cf81dSBarry Smith    Input Parameters:
1563e35cf81dSBarry Smith +  snes - the SNES context
1564e35cf81dSBarry Smith -  lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1565fe3ffe1eSBarry Smith          the Jacobian is built etc. -2 means rebuild at next chance but then never again
1566e35cf81dSBarry Smith 
1567e35cf81dSBarry Smith    Options Database Keys:
1568e35cf81dSBarry Smith .    -snes_lag_jacobian <lag>
1569e35cf81dSBarry Smith 
1570e35cf81dSBarry Smith    Notes:
1571e35cf81dSBarry Smith    The default is 1
1572e35cf81dSBarry Smith    The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1573fe3ffe1eSBarry Smith    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
1574fe3ffe1eSBarry Smith    at the next Newton step but never again (unless it is reset to another value)
1575e35cf81dSBarry Smith 
1576e35cf81dSBarry Smith    Level: intermediate
1577e35cf81dSBarry Smith 
1578e35cf81dSBarry Smith .keywords: SNES, nonlinear, set, convergence, tolerances
1579e35cf81dSBarry Smith 
1580e35cf81dSBarry Smith .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian()
1581e35cf81dSBarry Smith 
1582e35cf81dSBarry Smith @*/
15837087cfbeSBarry Smith PetscErrorCode  SNESSetLagJacobian(SNES snes,PetscInt lag)
1584e35cf81dSBarry Smith {
1585e35cf81dSBarry Smith   PetscFunctionBegin;
15860700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1587e32f2f54SBarry Smith   if (lag < -2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater");
1588e32f2f54SBarry Smith   if (!lag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1589c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(snes,lag,2);
1590e35cf81dSBarry Smith   snes->lagjacobian = lag;
1591e35cf81dSBarry Smith   PetscFunctionReturn(0);
1592e35cf81dSBarry Smith }
1593e35cf81dSBarry Smith 
1594e35cf81dSBarry Smith #undef __FUNCT__
1595e35cf81dSBarry Smith #define __FUNCT__ "SNESGetLagJacobian"
1596e35cf81dSBarry Smith /*@
1597e35cf81dSBarry Smith    SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt
1598e35cf81dSBarry Smith 
15993f9fe445SBarry Smith    Not Collective
1600e35cf81dSBarry Smith 
1601e35cf81dSBarry Smith    Input Parameter:
1602e35cf81dSBarry Smith .  snes - the SNES context
1603e35cf81dSBarry Smith 
1604e35cf81dSBarry Smith    Output Parameter:
1605e35cf81dSBarry Smith .   lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time
1606e35cf81dSBarry Smith          the Jacobian is built etc.
1607e35cf81dSBarry Smith 
1608e35cf81dSBarry Smith    Options Database Keys:
1609e35cf81dSBarry Smith .    -snes_lag_jacobian <lag>
1610e35cf81dSBarry Smith 
1611e35cf81dSBarry Smith    Notes:
1612e35cf81dSBarry Smith    The default is 1
1613e35cf81dSBarry Smith    The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1614e35cf81dSBarry Smith 
1615e35cf81dSBarry Smith    Level: intermediate
1616e35cf81dSBarry Smith 
1617e35cf81dSBarry Smith .keywords: SNES, nonlinear, set, convergence, tolerances
1618e35cf81dSBarry Smith 
1619e35cf81dSBarry Smith .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner()
1620e35cf81dSBarry Smith 
1621e35cf81dSBarry Smith @*/
16227087cfbeSBarry Smith PetscErrorCode  SNESGetLagJacobian(SNES snes,PetscInt *lag)
1623e35cf81dSBarry Smith {
1624e35cf81dSBarry Smith   PetscFunctionBegin;
16250700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1626e35cf81dSBarry Smith   *lag = snes->lagjacobian;
1627e35cf81dSBarry Smith   PetscFunctionReturn(0);
1628e35cf81dSBarry Smith }
1629e35cf81dSBarry Smith 
1630e35cf81dSBarry Smith #undef __FUNCT__
16314a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
16329b94acceSBarry Smith /*@
1633d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
16349b94acceSBarry Smith 
16353f9fe445SBarry Smith    Logically Collective on SNES
1636c7afd0dbSLois Curfman McInnes 
16379b94acceSBarry Smith    Input Parameters:
1638c7afd0dbSLois Curfman McInnes +  snes - the SNES context
163970441072SBarry Smith .  abstol - absolute convergence tolerance
164033174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
164133174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
164233174efeSLois Curfman McInnes            of the change in the solution between steps
164333174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1644c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1645fee21e36SBarry Smith 
164633174efeSLois Curfman McInnes    Options Database Keys:
164770441072SBarry Smith +    -snes_atol <abstol> - Sets abstol
1648c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1649c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1650c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1651c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
16529b94acceSBarry Smith 
1653d7a720efSLois Curfman McInnes    Notes:
16549b94acceSBarry Smith    The default maximum number of iterations is 50.
16559b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
16569b94acceSBarry Smith 
165736851e7fSLois Curfman McInnes    Level: intermediate
165836851e7fSLois Curfman McInnes 
165933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
16609b94acceSBarry Smith 
16612492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance()
16629b94acceSBarry Smith @*/
16637087cfbeSBarry Smith PetscErrorCode  SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
16649b94acceSBarry Smith {
16653a40ed3dSBarry Smith   PetscFunctionBegin;
16660700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1667c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,abstol,2);
1668c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,rtol,3);
1669c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,stol,4);
1670c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(snes,maxit,5);
1671c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(snes,maxf,6);
1672c5eb9154SBarry Smith 
167370441072SBarry Smith   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1674d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1675d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1676d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1677d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
16783a40ed3dSBarry Smith   PetscFunctionReturn(0);
16799b94acceSBarry Smith }
16809b94acceSBarry Smith 
16814a2ae208SSatish Balay #undef __FUNCT__
16824a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
16839b94acceSBarry Smith /*@
168433174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
168533174efeSLois Curfman McInnes 
1686c7afd0dbSLois Curfman McInnes    Not Collective
1687c7afd0dbSLois Curfman McInnes 
168833174efeSLois Curfman McInnes    Input Parameters:
1689c7afd0dbSLois Curfman McInnes +  snes - the SNES context
169085385478SLisandro Dalcin .  atol - absolute convergence tolerance
169133174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
169233174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
169333174efeSLois Curfman McInnes            of the change in the solution between steps
169433174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1695c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1696fee21e36SBarry Smith 
169733174efeSLois Curfman McInnes    Notes:
169833174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
169933174efeSLois Curfman McInnes 
170036851e7fSLois Curfman McInnes    Level: intermediate
170136851e7fSLois Curfman McInnes 
170233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
170333174efeSLois Curfman McInnes 
170433174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
170533174efeSLois Curfman McInnes @*/
17067087cfbeSBarry Smith PetscErrorCode  SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
170733174efeSLois Curfman McInnes {
17083a40ed3dSBarry Smith   PetscFunctionBegin;
17090700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
171085385478SLisandro Dalcin   if (atol)  *atol  = snes->abstol;
171133174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
171233174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
171333174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
171433174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
17153a40ed3dSBarry Smith   PetscFunctionReturn(0);
171633174efeSLois Curfman McInnes }
171733174efeSLois Curfman McInnes 
17184a2ae208SSatish Balay #undef __FUNCT__
17194a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
172033174efeSLois Curfman McInnes /*@
17219b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
17229b94acceSBarry Smith 
17233f9fe445SBarry Smith    Logically Collective on SNES
1724fee21e36SBarry Smith 
1725c7afd0dbSLois Curfman McInnes    Input Parameters:
1726c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1727c7afd0dbSLois Curfman McInnes -  tol - tolerance
1728c7afd0dbSLois Curfman McInnes 
17299b94acceSBarry Smith    Options Database Key:
1730c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
17319b94acceSBarry Smith 
173236851e7fSLois Curfman McInnes    Level: intermediate
173336851e7fSLois Curfman McInnes 
17349b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
17359b94acceSBarry Smith 
17362492ecdbSBarry Smith .seealso: SNESSetTolerances()
17379b94acceSBarry Smith @*/
17387087cfbeSBarry Smith PetscErrorCode  SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
17399b94acceSBarry Smith {
17403a40ed3dSBarry Smith   PetscFunctionBegin;
17410700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1742c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,tol,2);
17439b94acceSBarry Smith   snes->deltatol = tol;
17443a40ed3dSBarry Smith   PetscFunctionReturn(0);
17459b94acceSBarry Smith }
17469b94acceSBarry Smith 
1747df9fa365SBarry Smith /*
1748df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1749df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1750df9fa365SBarry Smith    macros instead of functions
1751df9fa365SBarry Smith */
17524a2ae208SSatish Balay #undef __FUNCT__
1753a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLG"
17547087cfbeSBarry Smith PetscErrorCode  SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1755ce1608b8SBarry Smith {
1756dfbe8321SBarry Smith   PetscErrorCode ierr;
1757ce1608b8SBarry Smith 
1758ce1608b8SBarry Smith   PetscFunctionBegin;
17590700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1760a6570f20SBarry Smith   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1761ce1608b8SBarry Smith   PetscFunctionReturn(0);
1762ce1608b8SBarry Smith }
1763ce1608b8SBarry Smith 
17644a2ae208SSatish Balay #undef __FUNCT__
1765a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLGCreate"
17667087cfbeSBarry Smith PetscErrorCode  SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1767df9fa365SBarry Smith {
1768dfbe8321SBarry Smith   PetscErrorCode ierr;
1769df9fa365SBarry Smith 
1770df9fa365SBarry Smith   PetscFunctionBegin;
1771a6570f20SBarry Smith   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1772df9fa365SBarry Smith   PetscFunctionReturn(0);
1773df9fa365SBarry Smith }
1774df9fa365SBarry Smith 
17754a2ae208SSatish Balay #undef __FUNCT__
1776a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLGDestroy"
17777087cfbeSBarry Smith PetscErrorCode  SNESMonitorLGDestroy(PetscDrawLG draw)
1778df9fa365SBarry Smith {
1779dfbe8321SBarry Smith   PetscErrorCode ierr;
1780df9fa365SBarry Smith 
1781df9fa365SBarry Smith   PetscFunctionBegin;
1782a6570f20SBarry Smith   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1783df9fa365SBarry Smith   PetscFunctionReturn(0);
1784df9fa365SBarry Smith }
1785df9fa365SBarry Smith 
17867087cfbeSBarry Smith extern PetscErrorCode  SNESMonitorRange_Private(SNES,PetscInt,PetscReal*);
1787b271bb04SBarry Smith #undef __FUNCT__
1788b271bb04SBarry Smith #define __FUNCT__ "SNESMonitorLGRange"
17897087cfbeSBarry Smith PetscErrorCode  SNESMonitorLGRange(SNES snes,PetscInt n,PetscReal rnorm,void *monctx)
1790b271bb04SBarry Smith {
1791b271bb04SBarry Smith   PetscDrawLG      lg;
1792b271bb04SBarry Smith   PetscErrorCode   ierr;
1793b271bb04SBarry Smith   PetscReal        x,y,per;
1794b271bb04SBarry Smith   PetscViewer      v = (PetscViewer)monctx;
1795b271bb04SBarry Smith   static PetscReal prev; /* should be in the context */
1796b271bb04SBarry Smith   PetscDraw        draw;
1797b271bb04SBarry Smith   PetscFunctionBegin;
1798b271bb04SBarry Smith   if (!monctx) {
1799b271bb04SBarry Smith     MPI_Comm    comm;
1800b271bb04SBarry Smith 
1801b271bb04SBarry Smith     ierr   = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
1802b271bb04SBarry Smith     v      = PETSC_VIEWER_DRAW_(comm);
1803b271bb04SBarry Smith   }
1804b271bb04SBarry Smith   ierr   = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
1805b271bb04SBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1806b271bb04SBarry Smith   ierr   = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1807b271bb04SBarry Smith   ierr   = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
1808b271bb04SBarry Smith   x = (PetscReal) n;
1809b271bb04SBarry Smith   if (rnorm > 0.0) y = log10(rnorm); else y = -15.0;
1810b271bb04SBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1811b271bb04SBarry Smith   if (n < 20 || !(n % 5)) {
1812b271bb04SBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1813b271bb04SBarry Smith   }
1814b271bb04SBarry Smith 
1815b271bb04SBarry Smith   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
1816b271bb04SBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1817b271bb04SBarry Smith   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1818b271bb04SBarry Smith   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
1819b271bb04SBarry Smith   ierr =  SNESMonitorRange_Private(snes,n,&per);CHKERRQ(ierr);
1820b271bb04SBarry Smith   x = (PetscReal) n;
1821b271bb04SBarry Smith   y = 100.0*per;
1822b271bb04SBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1823b271bb04SBarry Smith   if (n < 20 || !(n % 5)) {
1824b271bb04SBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1825b271bb04SBarry Smith   }
1826b271bb04SBarry Smith 
1827b271bb04SBarry Smith   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
1828b271bb04SBarry Smith   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1829b271bb04SBarry Smith   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1830b271bb04SBarry Smith   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
1831b271bb04SBarry Smith   x = (PetscReal) n;
1832b271bb04SBarry Smith   y = (prev - rnorm)/prev;
1833b271bb04SBarry Smith   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1834b271bb04SBarry Smith   if (n < 20 || !(n % 5)) {
1835b271bb04SBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1836b271bb04SBarry Smith   }
1837b271bb04SBarry Smith 
1838b271bb04SBarry Smith   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
1839b271bb04SBarry Smith   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
1840b271bb04SBarry Smith   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
1841b271bb04SBarry Smith   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
1842b271bb04SBarry Smith   x = (PetscReal) n;
1843b271bb04SBarry Smith   y = (prev - rnorm)/(prev*per);
1844b271bb04SBarry Smith   if (n > 2) { /*skip initial crazy value */
1845b271bb04SBarry Smith     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
1846b271bb04SBarry Smith   }
1847b271bb04SBarry Smith   if (n < 20 || !(n % 5)) {
1848b271bb04SBarry Smith     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
1849b271bb04SBarry Smith   }
1850b271bb04SBarry Smith   prev = rnorm;
1851b271bb04SBarry Smith   PetscFunctionReturn(0);
1852b271bb04SBarry Smith }
1853b271bb04SBarry Smith 
1854b271bb04SBarry Smith #undef __FUNCT__
1855b271bb04SBarry Smith #define __FUNCT__ "SNESMonitorLGRangeCreate"
18567087cfbeSBarry Smith PetscErrorCode  SNESMonitorLGRangeCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1857b271bb04SBarry Smith {
1858b271bb04SBarry Smith   PetscErrorCode ierr;
1859b271bb04SBarry Smith 
1860b271bb04SBarry Smith   PetscFunctionBegin;
1861b271bb04SBarry Smith   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1862b271bb04SBarry Smith   PetscFunctionReturn(0);
1863b271bb04SBarry Smith }
1864b271bb04SBarry Smith 
1865b271bb04SBarry Smith #undef __FUNCT__
1866b271bb04SBarry Smith #define __FUNCT__ "SNESMonitorLGRangeDestroy"
18677087cfbeSBarry Smith PetscErrorCode  SNESMonitorLGRangeDestroy(PetscDrawLG draw)
1868b271bb04SBarry Smith {
1869b271bb04SBarry Smith   PetscErrorCode ierr;
1870b271bb04SBarry Smith 
1871b271bb04SBarry Smith   PetscFunctionBegin;
1872b271bb04SBarry Smith   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1873b271bb04SBarry Smith   PetscFunctionReturn(0);
1874b271bb04SBarry Smith }
1875b271bb04SBarry Smith 
18769b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
18779b94acceSBarry Smith 
18784a2ae208SSatish Balay #undef __FUNCT__
1879a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorSet"
18809b94acceSBarry Smith /*@C
1881a6570f20SBarry Smith    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
18829b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
18839b94acceSBarry Smith    progress.
18849b94acceSBarry Smith 
18853f9fe445SBarry Smith    Logically Collective on SNES
1886fee21e36SBarry Smith 
1887c7afd0dbSLois Curfman McInnes    Input Parameters:
1888c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1889c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1890b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1891e8105e01SRichard Katz           monitor routine (use PETSC_NULL if no context is desired)
1892b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1893b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
18949b94acceSBarry Smith 
1895c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1896a7cc72afSBarry Smith $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1897c7afd0dbSLois Curfman McInnes 
1898c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1899c7afd0dbSLois Curfman McInnes .    its - iteration number
1900c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
190140a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
19029b94acceSBarry Smith 
19039665c990SLois Curfman McInnes    Options Database Keys:
1904a6570f20SBarry Smith +    -snes_monitor        - sets SNESMonitorDefault()
1905a6570f20SBarry Smith .    -snes_monitor_draw    - sets line graph monitor,
1906a6570f20SBarry Smith                             uses SNESMonitorLGCreate()
1907a6570f20SBarry Smith _    -snes_monitor_cancel - cancels all monitors that have
1908c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1909a6570f20SBarry Smith                             calls to SNESMonitorSet(), but
1910c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1911c7afd0dbSLois Curfman McInnes                             the options database.
19129665c990SLois Curfman McInnes 
1913639f9d9dSBarry Smith    Notes:
19146bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
1915a6570f20SBarry Smith    SNESMonitorSet() multiple times; all will be called in the
19166bc08f3fSLois Curfman McInnes    order in which they were set.
1917639f9d9dSBarry Smith 
1918025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each SNES object
1919025f1a04SBarry Smith 
192036851e7fSLois Curfman McInnes    Level: intermediate
192136851e7fSLois Curfman McInnes 
19229b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
19239b94acceSBarry Smith 
1924a6570f20SBarry Smith .seealso: SNESMonitorDefault(), SNESMonitorCancel()
19259b94acceSBarry Smith @*/
19267087cfbeSBarry Smith PetscErrorCode  SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
19279b94acceSBarry Smith {
1928b90d0a6eSBarry Smith   PetscInt i;
1929b90d0a6eSBarry Smith 
19303a40ed3dSBarry Smith   PetscFunctionBegin;
19310700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
193217186662SBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1933b90d0a6eSBarry Smith   for (i=0; i<snes->numbermonitors;i++) {
1934b90d0a6eSBarry Smith     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1935b90d0a6eSBarry Smith 
1936b90d0a6eSBarry Smith     /* check if both default monitors that share common ASCII viewer */
1937b90d0a6eSBarry Smith     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1938b90d0a6eSBarry Smith       if (mctx && snes->monitorcontext[i]) {
1939b90d0a6eSBarry Smith         PetscErrorCode          ierr;
1940b90d0a6eSBarry Smith         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1941b90d0a6eSBarry Smith         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1942b90d0a6eSBarry Smith         if (viewer1->viewer == viewer2->viewer) {
1943b90d0a6eSBarry Smith           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1944b90d0a6eSBarry Smith           PetscFunctionReturn(0);
1945b90d0a6eSBarry Smith         }
1946b90d0a6eSBarry Smith       }
1947b90d0a6eSBarry Smith     }
1948b90d0a6eSBarry Smith   }
1949b90d0a6eSBarry Smith   snes->monitor[snes->numbermonitors]           = monitor;
1950b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1951639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
19523a40ed3dSBarry Smith   PetscFunctionReturn(0);
19539b94acceSBarry Smith }
19549b94acceSBarry Smith 
19554a2ae208SSatish Balay #undef __FUNCT__
1956a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorCancel"
19575cd90555SBarry Smith /*@C
1958a6570f20SBarry Smith    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
19595cd90555SBarry Smith 
19603f9fe445SBarry Smith    Logically Collective on SNES
1961c7afd0dbSLois Curfman McInnes 
19625cd90555SBarry Smith    Input Parameters:
19635cd90555SBarry Smith .  snes - the SNES context
19645cd90555SBarry Smith 
19651a480d89SAdministrator    Options Database Key:
1966a6570f20SBarry Smith .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1967a6570f20SBarry Smith     into a code by calls to SNESMonitorSet(), but does not cancel those
1968c7afd0dbSLois Curfman McInnes     set via the options database
19695cd90555SBarry Smith 
19705cd90555SBarry Smith    Notes:
19715cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
19725cd90555SBarry Smith 
197336851e7fSLois Curfman McInnes    Level: intermediate
197436851e7fSLois Curfman McInnes 
19755cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
19765cd90555SBarry Smith 
1977a6570f20SBarry Smith .seealso: SNESMonitorDefault(), SNESMonitorSet()
19785cd90555SBarry Smith @*/
19797087cfbeSBarry Smith PetscErrorCode  SNESMonitorCancel(SNES snes)
19805cd90555SBarry Smith {
1981d952e501SBarry Smith   PetscErrorCode ierr;
1982d952e501SBarry Smith   PetscInt       i;
1983d952e501SBarry Smith 
19845cd90555SBarry Smith   PetscFunctionBegin;
19850700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1986d952e501SBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1987d952e501SBarry Smith     if (snes->monitordestroy[i]) {
1988d952e501SBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1989d952e501SBarry Smith     }
1990d952e501SBarry Smith   }
19915cd90555SBarry Smith   snes->numbermonitors = 0;
19925cd90555SBarry Smith   PetscFunctionReturn(0);
19935cd90555SBarry Smith }
19945cd90555SBarry Smith 
19954a2ae208SSatish Balay #undef __FUNCT__
19964a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
19979b94acceSBarry Smith /*@C
19989b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
19999b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
20009b94acceSBarry Smith 
20013f9fe445SBarry Smith    Logically Collective on SNES
2002fee21e36SBarry Smith 
2003c7afd0dbSLois Curfman McInnes    Input Parameters:
2004c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2005c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
20067f7931b9SBarry Smith .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
20077f7931b9SBarry Smith -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
20089b94acceSBarry Smith 
2009c7afd0dbSLois Curfman McInnes    Calling sequence of func:
201006ee9f85SBarry Smith $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
2011c7afd0dbSLois Curfman McInnes 
2012c7afd0dbSLois Curfman McInnes +    snes - the SNES context
201306ee9f85SBarry Smith .    it - current iteration (0 is the first and is before any Newton step)
2014c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
2015184914b5SBarry Smith .    reason - reason for convergence/divergence
2016c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
20174b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
20184b27c08aSLois Curfman McInnes -    f - 2-norm of function
20199b94acceSBarry Smith 
202036851e7fSLois Curfman McInnes    Level: advanced
202136851e7fSLois Curfman McInnes 
20229b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
20239b94acceSBarry Smith 
202485385478SLisandro Dalcin .seealso: SNESDefaultConverged(), SNESSkipConverged()
20259b94acceSBarry Smith @*/
20267087cfbeSBarry Smith PetscErrorCode  SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
20279b94acceSBarry Smith {
20287f7931b9SBarry Smith   PetscErrorCode ierr;
20297f7931b9SBarry Smith 
20303a40ed3dSBarry Smith   PetscFunctionBegin;
20310700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
203285385478SLisandro Dalcin   if (!func) func = SNESSkipConverged;
20337f7931b9SBarry Smith   if (snes->ops->convergeddestroy) {
20347f7931b9SBarry Smith     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
20357f7931b9SBarry Smith   }
203685385478SLisandro Dalcin   snes->ops->converged        = func;
20377f7931b9SBarry Smith   snes->ops->convergeddestroy = destroy;
203885385478SLisandro Dalcin   snes->cnvP                  = cctx;
20393a40ed3dSBarry Smith   PetscFunctionReturn(0);
20409b94acceSBarry Smith }
20419b94acceSBarry Smith 
20424a2ae208SSatish Balay #undef __FUNCT__
20434a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
204452baeb72SSatish Balay /*@
2045184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
2046184914b5SBarry Smith 
2047184914b5SBarry Smith    Not Collective
2048184914b5SBarry Smith 
2049184914b5SBarry Smith    Input Parameter:
2050184914b5SBarry Smith .  snes - the SNES context
2051184914b5SBarry Smith 
2052184914b5SBarry Smith    Output Parameter:
20534d0a8057SBarry Smith .  reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the
2054184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
2055184914b5SBarry Smith 
2056184914b5SBarry Smith    Level: intermediate
2057184914b5SBarry Smith 
2058184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
2059184914b5SBarry Smith 
2060184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
2061184914b5SBarry Smith 
206285385478SLisandro Dalcin .seealso: SNESSetConvergenceTest(), SNESConvergedReason
2063184914b5SBarry Smith @*/
20647087cfbeSBarry Smith PetscErrorCode  SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
2065184914b5SBarry Smith {
2066184914b5SBarry Smith   PetscFunctionBegin;
20670700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
20684482741eSBarry Smith   PetscValidPointer(reason,2);
2069184914b5SBarry Smith   *reason = snes->reason;
2070184914b5SBarry Smith   PetscFunctionReturn(0);
2071184914b5SBarry Smith }
2072184914b5SBarry Smith 
20734a2ae208SSatish Balay #undef __FUNCT__
20744a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
2075c9005455SLois Curfman McInnes /*@
2076c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
2077c9005455SLois Curfman McInnes 
20783f9fe445SBarry Smith    Logically Collective on SNES
2079fee21e36SBarry Smith 
2080c7afd0dbSLois Curfman McInnes    Input Parameters:
2081c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
20828c7482ecSBarry Smith .  a   - array to hold history, this array will contain the function norms computed at each step
2083cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
2084758f92a0SBarry Smith .  na  - size of a and its
208564731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
2086758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
2087c7afd0dbSLois Curfman McInnes 
2088c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
2089c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
2090c9005455SLois Curfman McInnes    during the section of code that is being timed.
2091c9005455SLois Curfman McInnes 
209236851e7fSLois Curfman McInnes    Level: intermediate
209336851e7fSLois Curfman McInnes 
2094c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
2095758f92a0SBarry Smith 
209608405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
2097758f92a0SBarry Smith 
2098c9005455SLois Curfman McInnes @*/
20997087cfbeSBarry Smith PetscErrorCode  SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscBool  reset)
2100c9005455SLois Curfman McInnes {
21013a40ed3dSBarry Smith   PetscFunctionBegin;
21020700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
21034482741eSBarry Smith   if (na)  PetscValidScalarPointer(a,2);
2104a562a398SLisandro Dalcin   if (its) PetscValidIntPointer(its,3);
2105c9005455SLois Curfman McInnes   snes->conv_hist       = a;
2106758f92a0SBarry Smith   snes->conv_hist_its   = its;
2107758f92a0SBarry Smith   snes->conv_hist_max   = na;
2108a12bc48fSLisandro Dalcin   snes->conv_hist_len   = 0;
2109758f92a0SBarry Smith   snes->conv_hist_reset = reset;
2110758f92a0SBarry Smith   PetscFunctionReturn(0);
2111758f92a0SBarry Smith }
2112758f92a0SBarry Smith 
21134a2ae208SSatish Balay #undef __FUNCT__
21144a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
21150c4c9dddSBarry Smith /*@C
2116758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
2117758f92a0SBarry Smith 
21183f9fe445SBarry Smith    Not Collective
2119758f92a0SBarry Smith 
2120758f92a0SBarry Smith    Input Parameter:
2121758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
2122758f92a0SBarry Smith 
2123758f92a0SBarry Smith    Output Parameters:
2124758f92a0SBarry Smith .  a   - array to hold history
2125758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
2126758f92a0SBarry Smith          negative if not converged) for each solve.
2127758f92a0SBarry Smith -  na  - size of a and its
2128758f92a0SBarry Smith 
2129758f92a0SBarry Smith    Notes:
2130758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
2131758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
2132758f92a0SBarry Smith 
2133758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
2134758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
2135758f92a0SBarry Smith    during the section of code that is being timed.
2136758f92a0SBarry Smith 
2137758f92a0SBarry Smith    Level: intermediate
2138758f92a0SBarry Smith 
2139758f92a0SBarry Smith .keywords: SNES, get, convergence, history
2140758f92a0SBarry Smith 
2141758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
2142758f92a0SBarry Smith 
2143758f92a0SBarry Smith @*/
21447087cfbeSBarry Smith PetscErrorCode  SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
2145758f92a0SBarry Smith {
2146758f92a0SBarry Smith   PetscFunctionBegin;
21470700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2148758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
2149758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
2150758f92a0SBarry Smith   if (na)  *na  = snes->conv_hist_len;
21513a40ed3dSBarry Smith   PetscFunctionReturn(0);
2152c9005455SLois Curfman McInnes }
2153c9005455SLois Curfman McInnes 
2154e74ef692SMatthew Knepley #undef __FUNCT__
2155e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
2156ac226902SBarry Smith /*@C
215776b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
2158eec8f646SJed Brown   at the beginning of every iteration of the nonlinear solve. Specifically
21597e4bb74cSBarry Smith   it is called just before the Jacobian is "evaluated".
216076b2cf59SMatthew Knepley 
21613f9fe445SBarry Smith   Logically Collective on SNES
216276b2cf59SMatthew Knepley 
216376b2cf59SMatthew Knepley   Input Parameters:
216476b2cf59SMatthew Knepley . snes - The nonlinear solver context
216576b2cf59SMatthew Knepley . func - The function
216676b2cf59SMatthew Knepley 
216776b2cf59SMatthew Knepley   Calling sequence of func:
2168b5d30489SBarry Smith . func (SNES snes, PetscInt step);
216976b2cf59SMatthew Knepley 
217076b2cf59SMatthew Knepley . step - The current step of the iteration
217176b2cf59SMatthew Knepley 
2172fe97e370SBarry Smith   Level: advanced
2173fe97e370SBarry Smith 
2174fe97e370SBarry Smith   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()
2175fe97e370SBarry Smith         This is not used by most users.
217676b2cf59SMatthew Knepley 
217776b2cf59SMatthew Knepley .keywords: SNES, update
2178b5d30489SBarry Smith 
217985385478SLisandro Dalcin .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
218076b2cf59SMatthew Knepley @*/
21817087cfbeSBarry Smith PetscErrorCode  SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
218276b2cf59SMatthew Knepley {
218376b2cf59SMatthew Knepley   PetscFunctionBegin;
21840700a824SBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID,1);
2185e7788613SBarry Smith   snes->ops->update = func;
218676b2cf59SMatthew Knepley   PetscFunctionReturn(0);
218776b2cf59SMatthew Knepley }
218876b2cf59SMatthew Knepley 
2189e74ef692SMatthew Knepley #undef __FUNCT__
2190e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
219176b2cf59SMatthew Knepley /*@
219276b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
219376b2cf59SMatthew Knepley 
219476b2cf59SMatthew Knepley   Not collective
219576b2cf59SMatthew Knepley 
219676b2cf59SMatthew Knepley   Input Parameters:
219776b2cf59SMatthew Knepley . snes - The nonlinear solver context
219876b2cf59SMatthew Knepley . step - The current step of the iteration
219976b2cf59SMatthew Knepley 
2200205452f4SMatthew Knepley   Level: intermediate
2201205452f4SMatthew Knepley 
220276b2cf59SMatthew Knepley .keywords: SNES, update
2203a6570f20SBarry Smith .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
220476b2cf59SMatthew Knepley @*/
22057087cfbeSBarry Smith PetscErrorCode  SNESDefaultUpdate(SNES snes, PetscInt step)
220676b2cf59SMatthew Knepley {
220776b2cf59SMatthew Knepley   PetscFunctionBegin;
220876b2cf59SMatthew Knepley   PetscFunctionReturn(0);
220976b2cf59SMatthew Knepley }
221076b2cf59SMatthew Knepley 
22114a2ae208SSatish Balay #undef __FUNCT__
22124a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
22139b94acceSBarry Smith /*
22149b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
22159b94acceSBarry Smith    positive parameter delta.
22169b94acceSBarry Smith 
22179b94acceSBarry Smith     Input Parameters:
2218c7afd0dbSLois Curfman McInnes +   snes - the SNES context
22199b94acceSBarry Smith .   y - approximate solution of linear system
22209b94acceSBarry Smith .   fnorm - 2-norm of current function
2221c7afd0dbSLois Curfman McInnes -   delta - trust region size
22229b94acceSBarry Smith 
22239b94acceSBarry Smith     Output Parameters:
2224c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
22259b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
22269b94acceSBarry Smith     region, and exceeds zero otherwise.
2227c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
22289b94acceSBarry Smith 
22299b94acceSBarry Smith     Note:
22304b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
22319b94acceSBarry Smith     is set to be the maximum allowable step size.
22329b94acceSBarry Smith 
22339b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
22349b94acceSBarry Smith */
2235dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
22369b94acceSBarry Smith {
2237064f8208SBarry Smith   PetscReal      nrm;
2238ea709b57SSatish Balay   PetscScalar    cnorm;
2239dfbe8321SBarry Smith   PetscErrorCode ierr;
22403a40ed3dSBarry Smith 
22413a40ed3dSBarry Smith   PetscFunctionBegin;
22420700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
22430700a824SBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,2);
2244c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
2245184914b5SBarry Smith 
2246064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
2247064f8208SBarry Smith   if (nrm > *delta) {
2248064f8208SBarry Smith      nrm = *delta/nrm;
2249064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
2250064f8208SBarry Smith      cnorm = nrm;
22512dcb1b2aSMatthew Knepley      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
22529b94acceSBarry Smith      *ynorm = *delta;
22539b94acceSBarry Smith   } else {
22549b94acceSBarry Smith      *gpnorm = 0.0;
2255064f8208SBarry Smith      *ynorm = nrm;
22569b94acceSBarry Smith   }
22573a40ed3dSBarry Smith   PetscFunctionReturn(0);
22589b94acceSBarry Smith }
22599b94acceSBarry Smith 
22604a2ae208SSatish Balay #undef __FUNCT__
22614a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
22626ce558aeSBarry Smith /*@C
2263f69a0ea3SMatthew Knepley    SNESSolve - Solves a nonlinear system F(x) = b.
2264f69a0ea3SMatthew Knepley    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
22659b94acceSBarry Smith 
2266c7afd0dbSLois Curfman McInnes    Collective on SNES
2267c7afd0dbSLois Curfman McInnes 
2268b2002411SLois Curfman McInnes    Input Parameters:
2269c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2270f69a0ea3SMatthew Knepley .  b - the constant part of the equation, or PETSC_NULL to use zero.
227185385478SLisandro Dalcin -  x - the solution vector.
22729b94acceSBarry Smith 
2273b2002411SLois Curfman McInnes    Notes:
22748ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
22758ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
22768ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
22778ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
22788ddd3da0SLois Curfman McInnes 
227936851e7fSLois Curfman McInnes    Level: beginner
228036851e7fSLois Curfman McInnes 
22819b94acceSBarry Smith .keywords: SNES, nonlinear, solve
22829b94acceSBarry Smith 
228385385478SLisandro Dalcin .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
22849b94acceSBarry Smith @*/
22857087cfbeSBarry Smith PetscErrorCode  SNESSolve(SNES snes,Vec b,Vec x)
22869b94acceSBarry Smith {
2287dfbe8321SBarry Smith   PetscErrorCode ierr;
2288ace3abfcSBarry Smith   PetscBool      flg;
2289eabae89aSBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
2290eabae89aSBarry Smith   PetscViewer    viewer;
2291052efed2SBarry Smith 
22923a40ed3dSBarry Smith   PetscFunctionBegin;
22930700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
22940700a824SBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
2295f69a0ea3SMatthew Knepley   PetscCheckSameComm(snes,1,x,3);
22960700a824SBarry Smith   if (b) PetscValidHeaderSpecific(b,VEC_CLASSID,2);
229785385478SLisandro Dalcin   if (b) PetscCheckSameComm(snes,1,b,2);
229885385478SLisandro Dalcin 
229985385478SLisandro Dalcin   /* set solution vector */
230085385478SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
230185385478SLisandro Dalcin   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
230285385478SLisandro Dalcin   snes->vec_sol = x;
230385385478SLisandro Dalcin   /* set afine vector if provided */
230485385478SLisandro Dalcin   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
230585385478SLisandro Dalcin   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
230685385478SLisandro Dalcin   snes->vec_rhs = b;
230785385478SLisandro Dalcin 
230885385478SLisandro Dalcin   if (!snes->vec_func && snes->vec_rhs) {
230985385478SLisandro Dalcin     ierr = VecDuplicate(b, &snes->vec_func);CHKERRQ(ierr);
231070e92668SMatthew Knepley   }
23113f149594SLisandro Dalcin 
231270e92668SMatthew Knepley   ierr = SNESSetUp(snes);CHKERRQ(ierr);
23133f149594SLisandro Dalcin 
2314abc0a331SBarry Smith   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
231550ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2316d5e45103SBarry Smith 
23173f149594SLisandro Dalcin   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
23184936397dSBarry Smith   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
231985385478SLisandro Dalcin   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
23204936397dSBarry Smith   if (snes->domainerror){
23214936397dSBarry Smith     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
23224936397dSBarry Smith     snes->domainerror = PETSC_FALSE;
23234936397dSBarry Smith   }
232417186662SBarry Smith   if (!snes->reason) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
23253f149594SLisandro Dalcin 
23267adad957SLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2327eabae89aSBarry Smith   if (flg && !PetscPreLoadingOn) {
23287adad957SLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2329eabae89aSBarry Smith     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2330eabae89aSBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2331eabae89aSBarry Smith   }
2332eabae89aSBarry Smith 
233390d69ab7SBarry Smith   flg  = PETSC_FALSE;
2334acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg,PETSC_NULL);CHKERRQ(ierr);
2335da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
23365968eb51SBarry Smith   if (snes->printreason) {
23375968eb51SBarry Smith     if (snes->reason > 0) {
23387adad957SLisandro Dalcin       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
23395968eb51SBarry Smith     } else {
23407adad957SLisandro Dalcin       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
23415968eb51SBarry Smith     }
23425968eb51SBarry Smith   }
23435968eb51SBarry Smith 
2344e113a28aSBarry Smith   if (snes->errorifnotconverged && snes->reason < 0) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged");
23453a40ed3dSBarry Smith   PetscFunctionReturn(0);
23469b94acceSBarry Smith }
23479b94acceSBarry Smith 
23489b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
23499b94acceSBarry Smith 
23504a2ae208SSatish Balay #undef __FUNCT__
23514a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
235282bf6240SBarry Smith /*@C
23534b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
23549b94acceSBarry Smith 
2355fee21e36SBarry Smith    Collective on SNES
2356fee21e36SBarry Smith 
2357c7afd0dbSLois Curfman McInnes    Input Parameters:
2358c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2359454a90a3SBarry Smith -  type - a known method
2360c7afd0dbSLois Curfman McInnes 
2361c7afd0dbSLois Curfman McInnes    Options Database Key:
2362454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
2363c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
2364ae12b187SLois Curfman McInnes 
23659b94acceSBarry Smith    Notes:
2366e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
23674b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
2368c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
23694b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
2370c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
23719b94acceSBarry Smith 
2372ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
2373ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
2374ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
2375ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
2376ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
2377ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
2378ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
2379ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
2380ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
2381b0a32e0cSBarry Smith   appropriate method.
238236851e7fSLois Curfman McInnes 
238336851e7fSLois Curfman McInnes   Level: intermediate
2384a703fe33SLois Curfman McInnes 
2385454a90a3SBarry Smith .keywords: SNES, set, type
2386435da068SBarry Smith 
2387435da068SBarry Smith .seealso: SNESType, SNESCreate()
2388435da068SBarry Smith 
23899b94acceSBarry Smith @*/
23907087cfbeSBarry Smith PetscErrorCode  SNESSetType(SNES snes,const SNESType type)
23919b94acceSBarry Smith {
2392dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
2393ace3abfcSBarry Smith   PetscBool      match;
23943a40ed3dSBarry Smith 
23953a40ed3dSBarry Smith   PetscFunctionBegin;
23960700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
23974482741eSBarry Smith   PetscValidCharPointer(type,2);
239882bf6240SBarry Smith 
23996831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
24000f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
240192ff6ae8SBarry Smith 
24027adad957SLisandro Dalcin   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2403e32f2f54SBarry Smith   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
240475396ef9SLisandro Dalcin   /* Destroy the previous private SNES context */
240575396ef9SLisandro Dalcin   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
240675396ef9SLisandro Dalcin   /* Reinitialize function pointers in SNESOps structure */
240775396ef9SLisandro Dalcin   snes->ops->setup          = 0;
240875396ef9SLisandro Dalcin   snes->ops->solve          = 0;
240975396ef9SLisandro Dalcin   snes->ops->view           = 0;
241075396ef9SLisandro Dalcin   snes->ops->setfromoptions = 0;
241175396ef9SLisandro Dalcin   snes->ops->destroy        = 0;
241275396ef9SLisandro Dalcin   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
241375396ef9SLisandro Dalcin   snes->setupcalled = PETSC_FALSE;
24143a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
2415454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
24169fb22e1aSBarry Smith #if defined(PETSC_HAVE_AMS)
24179fb22e1aSBarry Smith   if (PetscAMSPublishAll) {
24189fb22e1aSBarry Smith     ierr = PetscObjectAMSPublish((PetscObject)snes);CHKERRQ(ierr);
24199fb22e1aSBarry Smith   }
24209fb22e1aSBarry Smith #endif
24213a40ed3dSBarry Smith   PetscFunctionReturn(0);
24229b94acceSBarry Smith }
24239b94acceSBarry Smith 
2424a847f771SSatish Balay 
24259b94acceSBarry Smith /* --------------------------------------------------------------------- */
24264a2ae208SSatish Balay #undef __FUNCT__
24274a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
242852baeb72SSatish Balay /*@
24299b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2430f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
24319b94acceSBarry Smith 
2432fee21e36SBarry Smith    Not Collective
2433fee21e36SBarry Smith 
243436851e7fSLois Curfman McInnes    Level: advanced
243536851e7fSLois Curfman McInnes 
24369b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
24379b94acceSBarry Smith 
24389b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
24399b94acceSBarry Smith @*/
24407087cfbeSBarry Smith PetscErrorCode  SNESRegisterDestroy(void)
24419b94acceSBarry Smith {
2442dfbe8321SBarry Smith   PetscErrorCode ierr;
244382bf6240SBarry Smith 
24443a40ed3dSBarry Smith   PetscFunctionBegin;
24451441b1d3SBarry Smith   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
24464c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
24473a40ed3dSBarry Smith   PetscFunctionReturn(0);
24489b94acceSBarry Smith }
24499b94acceSBarry Smith 
24504a2ae208SSatish Balay #undef __FUNCT__
24514a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
24529b94acceSBarry Smith /*@C
24539a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
24549b94acceSBarry Smith 
2455c7afd0dbSLois Curfman McInnes    Not Collective
2456c7afd0dbSLois Curfman McInnes 
24579b94acceSBarry Smith    Input Parameter:
24584b0e389bSBarry Smith .  snes - nonlinear solver context
24599b94acceSBarry Smith 
24609b94acceSBarry Smith    Output Parameter:
24613a7fca6bSBarry Smith .  type - SNES method (a character string)
24629b94acceSBarry Smith 
246336851e7fSLois Curfman McInnes    Level: intermediate
246436851e7fSLois Curfman McInnes 
2465454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
24669b94acceSBarry Smith @*/
24677087cfbeSBarry Smith PetscErrorCode  SNESGetType(SNES snes,const SNESType *type)
24689b94acceSBarry Smith {
24693a40ed3dSBarry Smith   PetscFunctionBegin;
24700700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
24714482741eSBarry Smith   PetscValidPointer(type,2);
24727adad957SLisandro Dalcin   *type = ((PetscObject)snes)->type_name;
24733a40ed3dSBarry Smith   PetscFunctionReturn(0);
24749b94acceSBarry Smith }
24759b94acceSBarry Smith 
24764a2ae208SSatish Balay #undef __FUNCT__
24774a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
247852baeb72SSatish Balay /*@
24799b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
24809b94acceSBarry Smith    stored.
24819b94acceSBarry Smith 
2482c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2483c7afd0dbSLois Curfman McInnes 
24849b94acceSBarry Smith    Input Parameter:
24859b94acceSBarry Smith .  snes - the SNES context
24869b94acceSBarry Smith 
24879b94acceSBarry Smith    Output Parameter:
24889b94acceSBarry Smith .  x - the solution
24899b94acceSBarry Smith 
249070e92668SMatthew Knepley    Level: intermediate
249136851e7fSLois Curfman McInnes 
24929b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
24939b94acceSBarry Smith 
249485385478SLisandro Dalcin .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
24959b94acceSBarry Smith @*/
24967087cfbeSBarry Smith PetscErrorCode  SNESGetSolution(SNES snes,Vec *x)
24979b94acceSBarry Smith {
24983a40ed3dSBarry Smith   PetscFunctionBegin;
24990700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
25004482741eSBarry Smith   PetscValidPointer(x,2);
250185385478SLisandro Dalcin   *x = snes->vec_sol;
250270e92668SMatthew Knepley   PetscFunctionReturn(0);
250370e92668SMatthew Knepley }
250470e92668SMatthew Knepley 
250570e92668SMatthew Knepley #undef __FUNCT__
25064a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
250752baeb72SSatish Balay /*@
25089b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
25099b94acceSBarry Smith    stored.
25109b94acceSBarry Smith 
2511c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2512c7afd0dbSLois Curfman McInnes 
25139b94acceSBarry Smith    Input Parameter:
25149b94acceSBarry Smith .  snes - the SNES context
25159b94acceSBarry Smith 
25169b94acceSBarry Smith    Output Parameter:
25179b94acceSBarry Smith .  x - the solution update
25189b94acceSBarry Smith 
251936851e7fSLois Curfman McInnes    Level: advanced
252036851e7fSLois Curfman McInnes 
25219b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
25229b94acceSBarry Smith 
252385385478SLisandro Dalcin .seealso: SNESGetSolution(), SNESGetFunction()
25249b94acceSBarry Smith @*/
25257087cfbeSBarry Smith PetscErrorCode  SNESGetSolutionUpdate(SNES snes,Vec *x)
25269b94acceSBarry Smith {
25273a40ed3dSBarry Smith   PetscFunctionBegin;
25280700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
25294482741eSBarry Smith   PetscValidPointer(x,2);
253085385478SLisandro Dalcin   *x = snes->vec_sol_update;
25313a40ed3dSBarry Smith   PetscFunctionReturn(0);
25329b94acceSBarry Smith }
25339b94acceSBarry Smith 
25344a2ae208SSatish Balay #undef __FUNCT__
25354a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
25369b94acceSBarry Smith /*@C
25373638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
25389b94acceSBarry Smith 
2539c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2540c7afd0dbSLois Curfman McInnes 
25419b94acceSBarry Smith    Input Parameter:
25429b94acceSBarry Smith .  snes - the SNES context
25439b94acceSBarry Smith 
25449b94acceSBarry Smith    Output Parameter:
25457bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
254670e92668SMatthew Knepley .  func - the function (or PETSC_NULL)
254770e92668SMatthew Knepley -  ctx - the function context (or PETSC_NULL)
25489b94acceSBarry Smith 
254936851e7fSLois Curfman McInnes    Level: advanced
255036851e7fSLois Curfman McInnes 
2551a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
25529b94acceSBarry Smith 
25534b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
25549b94acceSBarry Smith @*/
25557087cfbeSBarry Smith PetscErrorCode  SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
25569b94acceSBarry Smith {
25573a40ed3dSBarry Smith   PetscFunctionBegin;
25580700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
255985385478SLisandro Dalcin   if (r)    *r    = snes->vec_func;
2560e7788613SBarry Smith   if (func) *func = snes->ops->computefunction;
256170e92668SMatthew Knepley   if (ctx)  *ctx  = snes->funP;
25623a40ed3dSBarry Smith   PetscFunctionReturn(0);
25639b94acceSBarry Smith }
25649b94acceSBarry Smith 
25654a2ae208SSatish Balay #undef __FUNCT__
25664a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
25673c7409f5SSatish Balay /*@C
25683c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2569d850072dSLois Curfman McInnes    SNES options in the database.
25703c7409f5SSatish Balay 
25713f9fe445SBarry Smith    Logically Collective on SNES
2572fee21e36SBarry Smith 
2573c7afd0dbSLois Curfman McInnes    Input Parameter:
2574c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2575c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2576c7afd0dbSLois Curfman McInnes 
2577d850072dSLois Curfman McInnes    Notes:
2578a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2579c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2580d850072dSLois Curfman McInnes 
258136851e7fSLois Curfman McInnes    Level: advanced
258236851e7fSLois Curfman McInnes 
25833c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2584a86d99e1SLois Curfman McInnes 
2585a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
25863c7409f5SSatish Balay @*/
25877087cfbeSBarry Smith PetscErrorCode  SNESSetOptionsPrefix(SNES snes,const char prefix[])
25883c7409f5SSatish Balay {
2589dfbe8321SBarry Smith   PetscErrorCode ierr;
25903c7409f5SSatish Balay 
25913a40ed3dSBarry Smith   PetscFunctionBegin;
25920700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2593639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
25941cee3971SBarry Smith   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
259594b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
25963a40ed3dSBarry Smith   PetscFunctionReturn(0);
25973c7409f5SSatish Balay }
25983c7409f5SSatish Balay 
25994a2ae208SSatish Balay #undef __FUNCT__
26004a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
26013c7409f5SSatish Balay /*@C
2602f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2603d850072dSLois Curfman McInnes    SNES options in the database.
26043c7409f5SSatish Balay 
26053f9fe445SBarry Smith    Logically Collective on SNES
2606fee21e36SBarry Smith 
2607c7afd0dbSLois Curfman McInnes    Input Parameters:
2608c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2609c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2610c7afd0dbSLois Curfman McInnes 
2611d850072dSLois Curfman McInnes    Notes:
2612a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2613c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2614d850072dSLois Curfman McInnes 
261536851e7fSLois Curfman McInnes    Level: advanced
261636851e7fSLois Curfman McInnes 
26173c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2618a86d99e1SLois Curfman McInnes 
2619a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
26203c7409f5SSatish Balay @*/
26217087cfbeSBarry Smith PetscErrorCode  SNESAppendOptionsPrefix(SNES snes,const char prefix[])
26223c7409f5SSatish Balay {
2623dfbe8321SBarry Smith   PetscErrorCode ierr;
26243c7409f5SSatish Balay 
26253a40ed3dSBarry Smith   PetscFunctionBegin;
26260700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2627639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
26281cee3971SBarry Smith   if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);}
262994b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
26303a40ed3dSBarry Smith   PetscFunctionReturn(0);
26313c7409f5SSatish Balay }
26323c7409f5SSatish Balay 
26334a2ae208SSatish Balay #undef __FUNCT__
26344a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
26359ab63eb5SSatish Balay /*@C
26363c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
26373c7409f5SSatish Balay    SNES options in the database.
26383c7409f5SSatish Balay 
2639c7afd0dbSLois Curfman McInnes    Not Collective
2640c7afd0dbSLois Curfman McInnes 
26413c7409f5SSatish Balay    Input Parameter:
26423c7409f5SSatish Balay .  snes - the SNES context
26433c7409f5SSatish Balay 
26443c7409f5SSatish Balay    Output Parameter:
26453c7409f5SSatish Balay .  prefix - pointer to the prefix string used
26463c7409f5SSatish Balay 
26474ef407dbSRichard Tran Mills    Notes: On the fortran side, the user should pass in a string 'prefix' of
26489ab63eb5SSatish Balay    sufficient length to hold the prefix.
26499ab63eb5SSatish Balay 
265036851e7fSLois Curfman McInnes    Level: advanced
265136851e7fSLois Curfman McInnes 
26523c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2653a86d99e1SLois Curfman McInnes 
2654a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
26553c7409f5SSatish Balay @*/
26567087cfbeSBarry Smith PetscErrorCode  SNESGetOptionsPrefix(SNES snes,const char *prefix[])
26573c7409f5SSatish Balay {
2658dfbe8321SBarry Smith   PetscErrorCode ierr;
26593c7409f5SSatish Balay 
26603a40ed3dSBarry Smith   PetscFunctionBegin;
26610700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2662639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
26633a40ed3dSBarry Smith   PetscFunctionReturn(0);
26643c7409f5SSatish Balay }
26653c7409f5SSatish Balay 
2666b2002411SLois Curfman McInnes 
26674a2ae208SSatish Balay #undef __FUNCT__
26684a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
26693cea93caSBarry Smith /*@C
26703cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
26713cea93caSBarry Smith 
26727f6c08e0SMatthew Knepley   Level: advanced
26733cea93caSBarry Smith @*/
26747087cfbeSBarry Smith PetscErrorCode  SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2675b2002411SLois Curfman McInnes {
2676e2d1d2b7SBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
2677dfbe8321SBarry Smith   PetscErrorCode ierr;
2678b2002411SLois Curfman McInnes 
2679b2002411SLois Curfman McInnes   PetscFunctionBegin;
2680b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2681c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2682b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2683b2002411SLois Curfman McInnes }
2684da9b6338SBarry Smith 
2685da9b6338SBarry Smith #undef __FUNCT__
2686da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
26877087cfbeSBarry Smith PetscErrorCode  SNESTestLocalMin(SNES snes)
2688da9b6338SBarry Smith {
2689dfbe8321SBarry Smith   PetscErrorCode ierr;
269077431f27SBarry Smith   PetscInt       N,i,j;
2691da9b6338SBarry Smith   Vec            u,uh,fh;
2692da9b6338SBarry Smith   PetscScalar    value;
2693da9b6338SBarry Smith   PetscReal      norm;
2694da9b6338SBarry Smith 
2695da9b6338SBarry Smith   PetscFunctionBegin;
2696da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2697da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2698da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2699da9b6338SBarry Smith 
2700da9b6338SBarry Smith   /* currently only works for sequential */
2701da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2702da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2703da9b6338SBarry Smith   for (i=0; i<N; i++) {
2704da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
270577431f27SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2706da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2707ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2708da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
27093ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2710da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
271177431f27SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2712da9b6338SBarry Smith       value = -value;
2713da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2714da9b6338SBarry Smith     }
2715da9b6338SBarry Smith   }
2716da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2717da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2718da9b6338SBarry Smith   PetscFunctionReturn(0);
2719da9b6338SBarry Smith }
272071f87433Sdalcinl 
272171f87433Sdalcinl #undef __FUNCT__
2722fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPSetUseEW"
272371f87433Sdalcinl /*@
2724fa9f3622SBarry Smith    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
272571f87433Sdalcinl    computing relative tolerance for linear solvers within an inexact
272671f87433Sdalcinl    Newton method.
272771f87433Sdalcinl 
27283f9fe445SBarry Smith    Logically Collective on SNES
272971f87433Sdalcinl 
273071f87433Sdalcinl    Input Parameters:
273171f87433Sdalcinl +  snes - SNES context
273271f87433Sdalcinl -  flag - PETSC_TRUE or PETSC_FALSE
273371f87433Sdalcinl 
273464ba62caSBarry Smith     Options Database:
273564ba62caSBarry Smith +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
273664ba62caSBarry Smith .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
273764ba62caSBarry Smith .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
273864ba62caSBarry Smith .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
273964ba62caSBarry Smith .  -snes_ksp_ew_gamma <gamma> - Sets gamma
274064ba62caSBarry Smith .  -snes_ksp_ew_alpha <alpha> - Sets alpha
274164ba62caSBarry Smith .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
274264ba62caSBarry Smith -  -snes_ksp_ew_threshold <threshold> - Sets threshold
274364ba62caSBarry Smith 
274471f87433Sdalcinl    Notes:
274571f87433Sdalcinl    Currently, the default is to use a constant relative tolerance for
274671f87433Sdalcinl    the inner linear solvers.  Alternatively, one can use the
274771f87433Sdalcinl    Eisenstat-Walker method, where the relative convergence tolerance
274871f87433Sdalcinl    is reset at each Newton iteration according progress of the nonlinear
274971f87433Sdalcinl    solver.
275071f87433Sdalcinl 
275171f87433Sdalcinl    Level: advanced
275271f87433Sdalcinl 
275371f87433Sdalcinl    Reference:
275471f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
275571f87433Sdalcinl    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
275671f87433Sdalcinl 
275771f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
275871f87433Sdalcinl 
2759fa9f3622SBarry Smith .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
276071f87433Sdalcinl @*/
27617087cfbeSBarry Smith PetscErrorCode  SNESKSPSetUseEW(SNES snes,PetscBool  flag)
276271f87433Sdalcinl {
276371f87433Sdalcinl   PetscFunctionBegin;
27640700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2765acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(snes,flag,2);
276671f87433Sdalcinl   snes->ksp_ewconv = flag;
276771f87433Sdalcinl   PetscFunctionReturn(0);
276871f87433Sdalcinl }
276971f87433Sdalcinl 
277071f87433Sdalcinl #undef __FUNCT__
2771fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPGetUseEW"
277271f87433Sdalcinl /*@
2773fa9f3622SBarry Smith    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
277471f87433Sdalcinl    for computing relative tolerance for linear solvers within an
277571f87433Sdalcinl    inexact Newton method.
277671f87433Sdalcinl 
277771f87433Sdalcinl    Not Collective
277871f87433Sdalcinl 
277971f87433Sdalcinl    Input Parameter:
278071f87433Sdalcinl .  snes - SNES context
278171f87433Sdalcinl 
278271f87433Sdalcinl    Output Parameter:
278371f87433Sdalcinl .  flag - PETSC_TRUE or PETSC_FALSE
278471f87433Sdalcinl 
278571f87433Sdalcinl    Notes:
278671f87433Sdalcinl    Currently, the default is to use a constant relative tolerance for
278771f87433Sdalcinl    the inner linear solvers.  Alternatively, one can use the
278871f87433Sdalcinl    Eisenstat-Walker method, where the relative convergence tolerance
278971f87433Sdalcinl    is reset at each Newton iteration according progress of the nonlinear
279071f87433Sdalcinl    solver.
279171f87433Sdalcinl 
279271f87433Sdalcinl    Level: advanced
279371f87433Sdalcinl 
279471f87433Sdalcinl    Reference:
279571f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
279671f87433Sdalcinl    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
279771f87433Sdalcinl 
279871f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
279971f87433Sdalcinl 
2800fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
280171f87433Sdalcinl @*/
28027087cfbeSBarry Smith PetscErrorCode  SNESKSPGetUseEW(SNES snes, PetscBool  *flag)
280371f87433Sdalcinl {
280471f87433Sdalcinl   PetscFunctionBegin;
28050700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
280671f87433Sdalcinl   PetscValidPointer(flag,2);
280771f87433Sdalcinl   *flag = snes->ksp_ewconv;
280871f87433Sdalcinl   PetscFunctionReturn(0);
280971f87433Sdalcinl }
281071f87433Sdalcinl 
281171f87433Sdalcinl #undef __FUNCT__
2812fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPSetParametersEW"
281371f87433Sdalcinl /*@
2814fa9f3622SBarry Smith    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
281571f87433Sdalcinl    convergence criteria for the linear solvers within an inexact
281671f87433Sdalcinl    Newton method.
281771f87433Sdalcinl 
28183f9fe445SBarry Smith    Logically Collective on SNES
281971f87433Sdalcinl 
282071f87433Sdalcinl    Input Parameters:
282171f87433Sdalcinl +    snes - SNES context
282271f87433Sdalcinl .    version - version 1, 2 (default is 2) or 3
282371f87433Sdalcinl .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
282471f87433Sdalcinl .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
282571f87433Sdalcinl .    gamma - multiplicative factor for version 2 rtol computation
282671f87433Sdalcinl              (0 <= gamma2 <= 1)
282771f87433Sdalcinl .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
282871f87433Sdalcinl .    alpha2 - power for safeguard
282971f87433Sdalcinl -    threshold - threshold for imposing safeguard (0 < threshold < 1)
283071f87433Sdalcinl 
283171f87433Sdalcinl    Note:
283271f87433Sdalcinl    Version 3 was contributed by Luis Chacon, June 2006.
283371f87433Sdalcinl 
283471f87433Sdalcinl    Use PETSC_DEFAULT to retain the default for any of the parameters.
283571f87433Sdalcinl 
283671f87433Sdalcinl    Level: advanced
283771f87433Sdalcinl 
283871f87433Sdalcinl    Reference:
283971f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
284071f87433Sdalcinl    inexact Newton method", Utah State University Math. Stat. Dept. Res.
284171f87433Sdalcinl    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
284271f87433Sdalcinl 
284371f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
284471f87433Sdalcinl 
2845fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
284671f87433Sdalcinl @*/
28477087cfbeSBarry Smith PetscErrorCode  SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
284871f87433Sdalcinl 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
284971f87433Sdalcinl {
2850fa9f3622SBarry Smith   SNESKSPEW *kctx;
285171f87433Sdalcinl   PetscFunctionBegin;
28520700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2853fa9f3622SBarry Smith   kctx = (SNESKSPEW*)snes->kspconvctx;
2854e32f2f54SBarry Smith   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
2855c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(snes,version,2);
2856c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,rtol_0,3);
2857c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,rtol_max,4);
2858c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,gamma,5);
2859c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,alpha,6);
2860c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,alpha2,7);
2861c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(snes,threshold,8);
286271f87433Sdalcinl 
286371f87433Sdalcinl   if (version != PETSC_DEFAULT)   kctx->version   = version;
286471f87433Sdalcinl   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
286571f87433Sdalcinl   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
286671f87433Sdalcinl   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
286771f87433Sdalcinl   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
286871f87433Sdalcinl   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
286971f87433Sdalcinl   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
287071f87433Sdalcinl 
287171f87433Sdalcinl   if (kctx->version < 1 || kctx->version > 3) {
2872e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
287371f87433Sdalcinl   }
287471f87433Sdalcinl   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
2875e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
287671f87433Sdalcinl   }
287771f87433Sdalcinl   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
2878e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
287971f87433Sdalcinl   }
288071f87433Sdalcinl   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
2881e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
288271f87433Sdalcinl   }
288371f87433Sdalcinl   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
2884e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
288571f87433Sdalcinl   }
288671f87433Sdalcinl   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
2887e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
288871f87433Sdalcinl   }
288971f87433Sdalcinl   PetscFunctionReturn(0);
289071f87433Sdalcinl }
289171f87433Sdalcinl 
289271f87433Sdalcinl #undef __FUNCT__
2893fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPGetParametersEW"
289471f87433Sdalcinl /*@
2895fa9f3622SBarry Smith    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
289671f87433Sdalcinl    convergence criteria for the linear solvers within an inexact
289771f87433Sdalcinl    Newton method.
289871f87433Sdalcinl 
289971f87433Sdalcinl    Not Collective
290071f87433Sdalcinl 
290171f87433Sdalcinl    Input Parameters:
290271f87433Sdalcinl      snes - SNES context
290371f87433Sdalcinl 
290471f87433Sdalcinl    Output Parameters:
290571f87433Sdalcinl +    version - version 1, 2 (default is 2) or 3
290671f87433Sdalcinl .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
290771f87433Sdalcinl .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
290871f87433Sdalcinl .    gamma - multiplicative factor for version 2 rtol computation
290971f87433Sdalcinl              (0 <= gamma2 <= 1)
291071f87433Sdalcinl .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
291171f87433Sdalcinl .    alpha2 - power for safeguard
291271f87433Sdalcinl -    threshold - threshold for imposing safeguard (0 < threshold < 1)
291371f87433Sdalcinl 
291471f87433Sdalcinl    Level: advanced
291571f87433Sdalcinl 
291671f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
291771f87433Sdalcinl 
2918fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
291971f87433Sdalcinl @*/
29207087cfbeSBarry Smith PetscErrorCode  SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
292171f87433Sdalcinl 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
292271f87433Sdalcinl {
2923fa9f3622SBarry Smith   SNESKSPEW *kctx;
292471f87433Sdalcinl   PetscFunctionBegin;
29250700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
2926fa9f3622SBarry Smith   kctx = (SNESKSPEW*)snes->kspconvctx;
2927e32f2f54SBarry Smith   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
292871f87433Sdalcinl   if(version)   *version   = kctx->version;
292971f87433Sdalcinl   if(rtol_0)    *rtol_0    = kctx->rtol_0;
293071f87433Sdalcinl   if(rtol_max)  *rtol_max  = kctx->rtol_max;
293171f87433Sdalcinl   if(gamma)     *gamma     = kctx->gamma;
293271f87433Sdalcinl   if(alpha)     *alpha     = kctx->alpha;
293371f87433Sdalcinl   if(alpha2)    *alpha2    = kctx->alpha2;
293471f87433Sdalcinl   if(threshold) *threshold = kctx->threshold;
293571f87433Sdalcinl   PetscFunctionReturn(0);
293671f87433Sdalcinl }
293771f87433Sdalcinl 
293871f87433Sdalcinl #undef __FUNCT__
2939fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPEW_PreSolve"
2940fa9f3622SBarry Smith static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
294171f87433Sdalcinl {
294271f87433Sdalcinl   PetscErrorCode ierr;
2943fa9f3622SBarry Smith   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
294471f87433Sdalcinl   PetscReal      rtol=PETSC_DEFAULT,stol;
294571f87433Sdalcinl 
294671f87433Sdalcinl   PetscFunctionBegin;
2947e32f2f54SBarry Smith   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
294871f87433Sdalcinl   if (!snes->iter) { /* first time in, so use the original user rtol */
294971f87433Sdalcinl     rtol = kctx->rtol_0;
295071f87433Sdalcinl   } else {
295171f87433Sdalcinl     if (kctx->version == 1) {
295271f87433Sdalcinl       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
295371f87433Sdalcinl       if (rtol < 0.0) rtol = -rtol;
295471f87433Sdalcinl       stol = pow(kctx->rtol_last,kctx->alpha2);
295571f87433Sdalcinl       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
295671f87433Sdalcinl     } else if (kctx->version == 2) {
295771f87433Sdalcinl       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
295871f87433Sdalcinl       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
295971f87433Sdalcinl       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
296071f87433Sdalcinl     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
296171f87433Sdalcinl       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
296271f87433Sdalcinl       /* safeguard: avoid sharp decrease of rtol */
296371f87433Sdalcinl       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
296471f87433Sdalcinl       stol = PetscMax(rtol,stol);
296571f87433Sdalcinl       rtol = PetscMin(kctx->rtol_0,stol);
296671f87433Sdalcinl       /* safeguard: avoid oversolving */
296771f87433Sdalcinl       stol = kctx->gamma*(snes->ttol)/snes->norm;
296871f87433Sdalcinl       stol = PetscMax(rtol,stol);
296971f87433Sdalcinl       rtol = PetscMin(kctx->rtol_0,stol);
2970e32f2f54SBarry Smith     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
297171f87433Sdalcinl   }
297271f87433Sdalcinl   /* safeguard: avoid rtol greater than one */
297371f87433Sdalcinl   rtol = PetscMin(rtol,kctx->rtol_max);
297471f87433Sdalcinl   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
297571f87433Sdalcinl   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
297671f87433Sdalcinl   PetscFunctionReturn(0);
297771f87433Sdalcinl }
297871f87433Sdalcinl 
297971f87433Sdalcinl #undef __FUNCT__
2980fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPEW_PostSolve"
2981fa9f3622SBarry Smith static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
298271f87433Sdalcinl {
298371f87433Sdalcinl   PetscErrorCode ierr;
2984fa9f3622SBarry Smith   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
298571f87433Sdalcinl   PCSide         pcside;
298671f87433Sdalcinl   Vec            lres;
298771f87433Sdalcinl 
298871f87433Sdalcinl   PetscFunctionBegin;
2989e32f2f54SBarry Smith   if (!kctx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
299071f87433Sdalcinl   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
299171f87433Sdalcinl   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
299271f87433Sdalcinl   if (kctx->version == 1) {
2993b037da10SBarry Smith     ierr = KSPGetPCSide(ksp,&pcside);CHKERRQ(ierr);
299471f87433Sdalcinl     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
299571f87433Sdalcinl       /* KSP residual is true linear residual */
299671f87433Sdalcinl       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
299771f87433Sdalcinl     } else {
299871f87433Sdalcinl       /* KSP residual is preconditioned residual */
299971f87433Sdalcinl       /* compute true linear residual norm */
300071f87433Sdalcinl       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
300171f87433Sdalcinl       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
300271f87433Sdalcinl       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
300371f87433Sdalcinl       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
300471f87433Sdalcinl       ierr = VecDestroy(lres);CHKERRQ(ierr);
300571f87433Sdalcinl     }
300671f87433Sdalcinl   }
300771f87433Sdalcinl   PetscFunctionReturn(0);
300871f87433Sdalcinl }
300971f87433Sdalcinl 
301071f87433Sdalcinl #undef __FUNCT__
301171f87433Sdalcinl #define __FUNCT__ "SNES_KSPSolve"
301271f87433Sdalcinl PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
301371f87433Sdalcinl {
301471f87433Sdalcinl   PetscErrorCode ierr;
301571f87433Sdalcinl 
301671f87433Sdalcinl   PetscFunctionBegin;
3017fa9f3622SBarry Smith   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
301871f87433Sdalcinl   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
3019fa9f3622SBarry Smith   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
302071f87433Sdalcinl   PetscFunctionReturn(0);
302171f87433Sdalcinl }
30226c699258SBarry Smith 
30236c699258SBarry Smith #undef __FUNCT__
30246c699258SBarry Smith #define __FUNCT__ "SNESSetDM"
30256c699258SBarry Smith /*@
30266c699258SBarry Smith    SNESSetDM - Sets the DM that may be used by some preconditioners
30276c699258SBarry Smith 
30283f9fe445SBarry Smith    Logically Collective on SNES
30296c699258SBarry Smith 
30306c699258SBarry Smith    Input Parameters:
30316c699258SBarry Smith +  snes - the preconditioner context
30326c699258SBarry Smith -  dm - the dm
30336c699258SBarry Smith 
30346c699258SBarry Smith    Level: intermediate
30356c699258SBarry Smith 
30366c699258SBarry Smith 
30376c699258SBarry Smith .seealso: SNESGetDM(), KSPSetDM(), KSPGetDM()
30386c699258SBarry Smith @*/
30397087cfbeSBarry Smith PetscErrorCode  SNESSetDM(SNES snes,DM dm)
30406c699258SBarry Smith {
30416c699258SBarry Smith   PetscErrorCode ierr;
3042345fed2cSBarry Smith   KSP            ksp;
30436c699258SBarry Smith 
30446c699258SBarry Smith   PetscFunctionBegin;
30450700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
304670663e4aSLisandro Dalcin   ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
30476c699258SBarry Smith   if (snes->dm) {ierr = DMDestroy(snes->dm);CHKERRQ(ierr);}
30486c699258SBarry Smith   snes->dm = dm;
3049345fed2cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
3050345fed2cSBarry Smith   ierr = KSPSetDM(ksp,dm);CHKERRQ(ierr);
3051f22f69f0SBarry Smith   ierr = KSPSetDMActive(ksp,PETSC_FALSE);CHKERRQ(ierr);
30526c699258SBarry Smith   PetscFunctionReturn(0);
30536c699258SBarry Smith }
30546c699258SBarry Smith 
30556c699258SBarry Smith #undef __FUNCT__
30566c699258SBarry Smith #define __FUNCT__ "SNESGetDM"
30576c699258SBarry Smith /*@
30586c699258SBarry Smith    SNESGetDM - Gets the DM that may be used by some preconditioners
30596c699258SBarry Smith 
30603f9fe445SBarry Smith    Not Collective but DM obtained is parallel on SNES
30616c699258SBarry Smith 
30626c699258SBarry Smith    Input Parameter:
30636c699258SBarry Smith . snes - the preconditioner context
30646c699258SBarry Smith 
30656c699258SBarry Smith    Output Parameter:
30666c699258SBarry Smith .  dm - the dm
30676c699258SBarry Smith 
30686c699258SBarry Smith    Level: intermediate
30696c699258SBarry Smith 
30706c699258SBarry Smith 
30716c699258SBarry Smith .seealso: SNESSetDM(), KSPSetDM(), KSPGetDM()
30726c699258SBarry Smith @*/
30737087cfbeSBarry Smith PetscErrorCode  SNESGetDM(SNES snes,DM *dm)
30746c699258SBarry Smith {
30756c699258SBarry Smith   PetscFunctionBegin;
30760700a824SBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
30776c699258SBarry Smith   *dm = snes->dm;
30786c699258SBarry Smith   PetscFunctionReturn(0);
30796c699258SBarry Smith }
30800807856dSBarry Smith 
308169b4f73cSBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
3082e650e774SBarry Smith #include "mex.h"
308369b4f73cSBarry Smith 
30848f6e6473SBarry Smith typedef struct {char *funcname; mxArray *ctx;} SNESMatlabContext;
30858f6e6473SBarry Smith 
30860807856dSBarry Smith #undef __FUNCT__
30870807856dSBarry Smith #define __FUNCT__ "SNESComputeFunction_Matlab"
30880807856dSBarry Smith /*
30890807856dSBarry Smith    SNESComputeFunction_Matlab - Calls the function that has been set with
30900807856dSBarry Smith                          SNESSetFunctionMatlab().
30910807856dSBarry Smith 
30920807856dSBarry Smith    Collective on SNES
30930807856dSBarry Smith 
30940807856dSBarry Smith    Input Parameters:
30950807856dSBarry Smith +  snes - the SNES context
30960807856dSBarry Smith -  x - input vector
30970807856dSBarry Smith 
30980807856dSBarry Smith    Output Parameter:
30990807856dSBarry Smith .  y - function vector, as set by SNESSetFunction()
31000807856dSBarry Smith 
31010807856dSBarry Smith    Notes:
31020807856dSBarry Smith    SNESComputeFunction() is typically used within nonlinear solvers
31030807856dSBarry Smith    implementations, so most users would not generally call this routine
31040807856dSBarry Smith    themselves.
31050807856dSBarry Smith 
31060807856dSBarry Smith    Level: developer
31070807856dSBarry Smith 
31080807856dSBarry Smith .keywords: SNES, nonlinear, compute, function
31090807856dSBarry Smith 
31100807856dSBarry Smith .seealso: SNESSetFunction(), SNESGetFunction()
311161b2408cSBarry Smith */
31127087cfbeSBarry Smith PetscErrorCode  SNESComputeFunction_Matlab(SNES snes,Vec x,Vec y, void *ctx)
31130807856dSBarry Smith {
3114e650e774SBarry Smith   PetscErrorCode    ierr;
31158f6e6473SBarry Smith   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
31168f6e6473SBarry Smith   int               nlhs = 1,nrhs = 5;
31178f6e6473SBarry Smith   mxArray	    *plhs[1],*prhs[5];
311891621f2eSBarry Smith   long long int     lx = 0,ly = 0,ls = 0;
3119e650e774SBarry Smith 
31200807856dSBarry Smith   PetscFunctionBegin;
31210807856dSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
31220807856dSBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
31230807856dSBarry Smith   PetscValidHeaderSpecific(y,VEC_CLASSID,3);
31240807856dSBarry Smith   PetscCheckSameComm(snes,1,x,2);
31250807856dSBarry Smith   PetscCheckSameComm(snes,1,y,3);
31260807856dSBarry Smith 
31270807856dSBarry Smith   /* call Matlab function in ctx with arguments x and y */
3128e650e774SBarry Smith 
312991621f2eSBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3130e650e774SBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3131e650e774SBarry Smith   ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr);
313291621f2eSBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
313391621f2eSBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
313491621f2eSBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)ly);
31358f6e6473SBarry Smith   prhs[3] =  mxCreateString(sctx->funcname);
31368f6e6473SBarry Smith   prhs[4] =  sctx->ctx;
3137b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeFunctionInternal");CHKERRQ(ierr);
3138e650e774SBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3139e650e774SBarry Smith   mxDestroyArray(prhs[0]);
3140e650e774SBarry Smith   mxDestroyArray(prhs[1]);
3141e650e774SBarry Smith   mxDestroyArray(prhs[2]);
31428f6e6473SBarry Smith   mxDestroyArray(prhs[3]);
3143e650e774SBarry Smith   mxDestroyArray(plhs[0]);
31440807856dSBarry Smith   PetscFunctionReturn(0);
31450807856dSBarry Smith }
31460807856dSBarry Smith 
31470807856dSBarry Smith 
31480807856dSBarry Smith #undef __FUNCT__
31490807856dSBarry Smith #define __FUNCT__ "SNESSetFunctionMatlab"
315061b2408cSBarry Smith /*
31510807856dSBarry Smith    SNESSetFunctionMatlab - Sets the function evaluation routine and function
31520807856dSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
31530807856dSBarry Smith    equations from Matlab. Here the function is a string containing the name of a Matlab function
31540807856dSBarry Smith 
31550807856dSBarry Smith    Logically Collective on SNES
31560807856dSBarry Smith 
31570807856dSBarry Smith    Input Parameters:
31580807856dSBarry Smith +  snes - the SNES context
31590807856dSBarry Smith .  r - vector to store function value
31600807856dSBarry Smith -  func - function evaluation routine
31610807856dSBarry Smith 
31620807856dSBarry Smith    Calling sequence of func:
316361b2408cSBarry Smith $    func (SNES snes,Vec x,Vec f,void *ctx);
31640807856dSBarry Smith 
31650807856dSBarry Smith 
31660807856dSBarry Smith    Notes:
31670807856dSBarry Smith    The Newton-like methods typically solve linear systems of the form
31680807856dSBarry Smith $      f'(x) x = -f(x),
31690807856dSBarry Smith    where f'(x) denotes the Jacobian matrix and f(x) is the function.
31700807856dSBarry Smith 
31710807856dSBarry Smith    Level: beginner
31720807856dSBarry Smith 
31730807856dSBarry Smith .keywords: SNES, nonlinear, set, function
31740807856dSBarry Smith 
31750807856dSBarry Smith .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
317661b2408cSBarry Smith */
31777087cfbeSBarry Smith PetscErrorCode  SNESSetFunctionMatlab(SNES snes,Vec r,const char *func,mxArray *ctx)
31780807856dSBarry Smith {
31790807856dSBarry Smith   PetscErrorCode    ierr;
31808f6e6473SBarry Smith   SNESMatlabContext *sctx;
31810807856dSBarry Smith 
31820807856dSBarry Smith   PetscFunctionBegin;
31838f6e6473SBarry Smith   /* currently sctx is memory bleed */
31848f6e6473SBarry Smith   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
31858f6e6473SBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
31868f6e6473SBarry Smith   /*
31878f6e6473SBarry Smith      This should work, but it doesn't
31888f6e6473SBarry Smith   sctx->ctx = ctx;
31898f6e6473SBarry Smith   mexMakeArrayPersistent(sctx->ctx);
31908f6e6473SBarry Smith   */
31918f6e6473SBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
31928f6e6473SBarry Smith   ierr = SNESSetFunction(snes,r,SNESComputeFunction_Matlab,sctx);CHKERRQ(ierr);
31930807856dSBarry Smith   PetscFunctionReturn(0);
31940807856dSBarry Smith }
319569b4f73cSBarry Smith 
319661b2408cSBarry Smith #undef __FUNCT__
319761b2408cSBarry Smith #define __FUNCT__ "SNESComputeJacobian_Matlab"
319861b2408cSBarry Smith /*
319961b2408cSBarry Smith    SNESComputeJacobian_Matlab - Calls the function that has been set with
320061b2408cSBarry Smith                          SNESSetJacobianMatlab().
320161b2408cSBarry Smith 
320261b2408cSBarry Smith    Collective on SNES
320361b2408cSBarry Smith 
320461b2408cSBarry Smith    Input Parameters:
320561b2408cSBarry Smith +  snes - the SNES context
320661b2408cSBarry Smith .  x - input vector
320761b2408cSBarry Smith .  A, B - the matrices
320861b2408cSBarry Smith -  ctx - user context
320961b2408cSBarry Smith 
321061b2408cSBarry Smith    Output Parameter:
321161b2408cSBarry Smith .  flag - structure of the matrix
321261b2408cSBarry Smith 
321361b2408cSBarry Smith    Level: developer
321461b2408cSBarry Smith 
321561b2408cSBarry Smith .keywords: SNES, nonlinear, compute, function
321661b2408cSBarry Smith 
321761b2408cSBarry Smith .seealso: SNESSetFunction(), SNESGetFunction()
321861b2408cSBarry Smith @*/
32197087cfbeSBarry Smith PetscErrorCode  SNESComputeJacobian_Matlab(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag, void *ctx)
322061b2408cSBarry Smith {
322161b2408cSBarry Smith   PetscErrorCode    ierr;
322261b2408cSBarry Smith   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
322361b2408cSBarry Smith   int               nlhs = 2,nrhs = 6;
322461b2408cSBarry Smith   mxArray	    *plhs[2],*prhs[6];
322561b2408cSBarry Smith   long long int     lx = 0,lA = 0,ls = 0, lB = 0;
322661b2408cSBarry Smith 
322761b2408cSBarry Smith   PetscFunctionBegin;
322861b2408cSBarry Smith   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
322961b2408cSBarry Smith   PetscValidHeaderSpecific(x,VEC_CLASSID,2);
323061b2408cSBarry Smith 
323161b2408cSBarry Smith   /* call Matlab function in ctx with arguments x and y */
323261b2408cSBarry Smith 
323361b2408cSBarry Smith   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
323461b2408cSBarry Smith   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
323561b2408cSBarry Smith   ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr);
323661b2408cSBarry Smith   ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr);
323761b2408cSBarry Smith   prhs[0] =  mxCreateDoubleScalar((double)ls);
323861b2408cSBarry Smith   prhs[1] =  mxCreateDoubleScalar((double)lx);
323961b2408cSBarry Smith   prhs[2] =  mxCreateDoubleScalar((double)lA);
324061b2408cSBarry Smith   prhs[3] =  mxCreateDoubleScalar((double)lB);
324161b2408cSBarry Smith   prhs[4] =  mxCreateString(sctx->funcname);
324261b2408cSBarry Smith   prhs[5] =  sctx->ctx;
3243b807a863SBarry Smith   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESComputeJacobianInternal");CHKERRQ(ierr);
324461b2408cSBarry Smith   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
324561b2408cSBarry Smith   *flag   =  (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr);
324661b2408cSBarry Smith   mxDestroyArray(prhs[0]);
324761b2408cSBarry Smith   mxDestroyArray(prhs[1]);
324861b2408cSBarry Smith   mxDestroyArray(prhs[2]);
324961b2408cSBarry Smith   mxDestroyArray(prhs[3]);
325061b2408cSBarry Smith   mxDestroyArray(prhs[4]);
325161b2408cSBarry Smith   mxDestroyArray(plhs[0]);
325261b2408cSBarry Smith   mxDestroyArray(plhs[1]);
325361b2408cSBarry Smith   PetscFunctionReturn(0);
325461b2408cSBarry Smith }
325561b2408cSBarry Smith 
325661b2408cSBarry Smith 
325761b2408cSBarry Smith #undef __FUNCT__
325861b2408cSBarry Smith #define __FUNCT__ "SNESSetJacobianMatlab"
325961b2408cSBarry Smith /*
326061b2408cSBarry Smith    SNESSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices
326161b2408cSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
326261b2408cSBarry Smith    equations from Matlab. Here the function is a string containing the name of a Matlab function
326361b2408cSBarry Smith 
326461b2408cSBarry Smith    Logically Collective on SNES
326561b2408cSBarry Smith 
326661b2408cSBarry Smith    Input Parameters:
326761b2408cSBarry Smith +  snes - the SNES context
326861b2408cSBarry Smith .  A,B - Jacobian matrices
326961b2408cSBarry Smith .  func - function evaluation routine
327061b2408cSBarry Smith -  ctx - user context
327161b2408cSBarry Smith 
327261b2408cSBarry Smith    Calling sequence of func:
327361b2408cSBarry Smith $    flag = func (SNES snes,Vec x,Mat A,Mat B,void *ctx);
327461b2408cSBarry Smith 
327561b2408cSBarry Smith 
327661b2408cSBarry Smith    Level: developer
327761b2408cSBarry Smith 
327861b2408cSBarry Smith .keywords: SNES, nonlinear, set, function
327961b2408cSBarry Smith 
328061b2408cSBarry Smith .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
328161b2408cSBarry Smith */
32827087cfbeSBarry Smith PetscErrorCode  SNESSetJacobianMatlab(SNES snes,Mat A,Mat B,const char *func,mxArray *ctx)
328361b2408cSBarry Smith {
328461b2408cSBarry Smith   PetscErrorCode    ierr;
328561b2408cSBarry Smith   SNESMatlabContext *sctx;
328661b2408cSBarry Smith 
328761b2408cSBarry Smith   PetscFunctionBegin;
328861b2408cSBarry Smith   /* currently sctx is memory bleed */
328961b2408cSBarry Smith   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
329061b2408cSBarry Smith   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
329161b2408cSBarry Smith   /*
329261b2408cSBarry Smith      This should work, but it doesn't
329361b2408cSBarry Smith   sctx->ctx = ctx;
329461b2408cSBarry Smith   mexMakeArrayPersistent(sctx->ctx);
329561b2408cSBarry Smith   */
329661b2408cSBarry Smith   sctx->ctx = mxDuplicateArray(ctx);
329761b2408cSBarry Smith   ierr = SNESSetJacobian(snes,A,B,SNESComputeJacobian_Matlab,sctx);CHKERRQ(ierr);
329861b2408cSBarry Smith   PetscFunctionReturn(0);
329961b2408cSBarry Smith }
330069b4f73cSBarry Smith 
3301f9eb7ae2SShri Abhyankar #undef __FUNCT__
3302f9eb7ae2SShri Abhyankar #define __FUNCT__ "SNESMonitor_Matlab"
3303f9eb7ae2SShri Abhyankar /*
3304f9eb7ae2SShri Abhyankar    SNESMonitor_Matlab - Calls the function that has been set with SNESMonitorSetMatlab().
3305f9eb7ae2SShri Abhyankar 
3306f9eb7ae2SShri Abhyankar    Collective on SNES
3307f9eb7ae2SShri Abhyankar 
3308f9eb7ae2SShri Abhyankar .seealso: SNESSetFunction(), SNESGetFunction()
3309f9eb7ae2SShri Abhyankar @*/
33107087cfbeSBarry Smith PetscErrorCode  SNESMonitor_Matlab(SNES snes,PetscInt it, PetscReal fnorm, void *ctx)
3311f9eb7ae2SShri Abhyankar {
3312f9eb7ae2SShri Abhyankar   PetscErrorCode  ierr;
331348f37e70SShri Abhyankar   SNESMatlabContext *sctx = (SNESMatlabContext *)ctx;
3314f9eb7ae2SShri Abhyankar   int             nlhs = 1,nrhs = 6;
3315f9eb7ae2SShri Abhyankar   mxArray	  *plhs[1],*prhs[6];
3316f9eb7ae2SShri Abhyankar   long long int   lx = 0,ls = 0;
3317f9eb7ae2SShri Abhyankar   Vec             x=snes->vec_sol;
3318f9eb7ae2SShri Abhyankar 
3319f9eb7ae2SShri Abhyankar   PetscFunctionBegin;
3320f9eb7ae2SShri Abhyankar   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
3321f9eb7ae2SShri Abhyankar 
3322f9eb7ae2SShri Abhyankar   ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr);
3323f9eb7ae2SShri Abhyankar   ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr);
3324f9eb7ae2SShri Abhyankar   prhs[0] =  mxCreateDoubleScalar((double)ls);
3325f9eb7ae2SShri Abhyankar   prhs[1] =  mxCreateDoubleScalar((double)it);
3326f9eb7ae2SShri Abhyankar   prhs[2] =  mxCreateDoubleScalar((double)fnorm);
3327f9eb7ae2SShri Abhyankar   prhs[3] =  mxCreateDoubleScalar((double)lx);
3328f9eb7ae2SShri Abhyankar   prhs[4] =  mxCreateString(sctx->funcname);
3329f9eb7ae2SShri Abhyankar   prhs[5] =  sctx->ctx;
3330f9eb7ae2SShri Abhyankar   ierr    =  mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscSNESMonitorInternal");CHKERRQ(ierr);
3331f9eb7ae2SShri Abhyankar   ierr    =  mxGetScalar(plhs[0]);CHKERRQ(ierr);
3332f9eb7ae2SShri Abhyankar   mxDestroyArray(prhs[0]);
3333f9eb7ae2SShri Abhyankar   mxDestroyArray(prhs[1]);
3334f9eb7ae2SShri Abhyankar   mxDestroyArray(prhs[2]);
3335f9eb7ae2SShri Abhyankar   mxDestroyArray(prhs[3]);
3336f9eb7ae2SShri Abhyankar   mxDestroyArray(prhs[4]);
3337f9eb7ae2SShri Abhyankar   mxDestroyArray(plhs[0]);
3338f9eb7ae2SShri Abhyankar   PetscFunctionReturn(0);
3339f9eb7ae2SShri Abhyankar }
3340f9eb7ae2SShri Abhyankar 
3341f9eb7ae2SShri Abhyankar 
3342f9eb7ae2SShri Abhyankar #undef __FUNCT__
3343f9eb7ae2SShri Abhyankar #define __FUNCT__ "SNESMonitorSetMatlab"
3344f9eb7ae2SShri Abhyankar /*
3345f9eb7ae2SShri Abhyankar    SNESMonitorSetMatlab - Sets the monitor function from Matlab
3346f9eb7ae2SShri Abhyankar 
3347f9eb7ae2SShri Abhyankar    Level: developer
3348f9eb7ae2SShri Abhyankar 
3349f9eb7ae2SShri Abhyankar .keywords: SNES, nonlinear, set, function
3350f9eb7ae2SShri Abhyankar 
3351f9eb7ae2SShri Abhyankar .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
3352f9eb7ae2SShri Abhyankar */
33537087cfbeSBarry Smith PetscErrorCode  SNESMonitorSetMatlab(SNES snes,const char *func,mxArray *ctx)
3354f9eb7ae2SShri Abhyankar {
3355f9eb7ae2SShri Abhyankar   PetscErrorCode    ierr;
3356f9eb7ae2SShri Abhyankar   SNESMatlabContext *sctx;
3357f9eb7ae2SShri Abhyankar 
3358f9eb7ae2SShri Abhyankar   PetscFunctionBegin;
3359f9eb7ae2SShri Abhyankar   /* currently sctx is memory bleed */
3360f9eb7ae2SShri Abhyankar   ierr = PetscMalloc(sizeof(SNESMatlabContext),&sctx);CHKERRQ(ierr);
3361f9eb7ae2SShri Abhyankar   ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr);
3362f9eb7ae2SShri Abhyankar   /*
3363f9eb7ae2SShri Abhyankar      This should work, but it doesn't
3364f9eb7ae2SShri Abhyankar   sctx->ctx = ctx;
3365f9eb7ae2SShri Abhyankar   mexMakeArrayPersistent(sctx->ctx);
3366f9eb7ae2SShri Abhyankar   */
3367f9eb7ae2SShri Abhyankar   sctx->ctx = mxDuplicateArray(ctx);
3368f9eb7ae2SShri Abhyankar   ierr = SNESMonitorSet(snes,SNESMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr);
3369f9eb7ae2SShri Abhyankar   PetscFunctionReturn(0);
3370f9eb7ae2SShri Abhyankar }
3371f9eb7ae2SShri Abhyankar 
337269b4f73cSBarry Smith #endif
3373