xref: /petsc/src/snes/interface/snes.c (revision a805402700a16671df1c4d2ca73e243b5e674cf2)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
29b94acceSBarry Smith 
3b9147fbbSdalcinl #include "include/private/snesimpl.h"      /*I "petscsnes.h"  I*/
49b94acceSBarry Smith 
54c49b128SBarry Smith PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
68ba1e511SMatthew Knepley PetscFList SNESList              = PETSC_NULL;
78ba1e511SMatthew Knepley 
88ba1e511SMatthew Knepley /* Logging support */
9166c7f25SBarry Smith PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE;
10166c7f25SBarry Smith PetscLogEvent  SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval;
11a09944afSBarry Smith 
12a09944afSBarry Smith #undef __FUNCT__
134936397dSBarry Smith #define __FUNCT__ "SNESSetFunctionDomainError"
14e725d27bSBarry Smith /*@
154936397dSBarry Smith    SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not
164936397dSBarry Smith      in the functions domain. For example, negative pressure.
174936397dSBarry Smith 
184936397dSBarry Smith    Collective on SNES
194936397dSBarry Smith 
204936397dSBarry Smith    Input Parameters:
214936397dSBarry Smith .  SNES - the SNES context
224936397dSBarry Smith 
2328529972SSatish Balay    Level: advanced
244936397dSBarry Smith 
254936397dSBarry Smith .keywords: SNES, view
264936397dSBarry Smith 
274936397dSBarry Smith .seealso: SNESCreate(), SNESSetFunction()
284936397dSBarry Smith @*/
294936397dSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionDomainError(SNES snes)
304936397dSBarry Smith {
314936397dSBarry Smith   PetscFunctionBegin;
324936397dSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
334936397dSBarry Smith   snes->domainerror = PETSC_TRUE;
344936397dSBarry Smith   PetscFunctionReturn(0);
354936397dSBarry Smith }
364936397dSBarry Smith 
374936397dSBarry Smith #undef __FUNCT__
384a2ae208SSatish Balay #define __FUNCT__ "SNESView"
397e2c5f70SBarry Smith /*@C
409b94acceSBarry Smith    SNESView - Prints the SNES data structure.
419b94acceSBarry Smith 
424c49b128SBarry Smith    Collective on SNES
43fee21e36SBarry Smith 
44c7afd0dbSLois Curfman McInnes    Input Parameters:
45c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
46c7afd0dbSLois Curfman McInnes -  viewer - visualization context
47c7afd0dbSLois Curfman McInnes 
489b94acceSBarry Smith    Options Database Key:
49c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
509b94acceSBarry Smith 
519b94acceSBarry Smith    Notes:
529b94acceSBarry Smith    The available visualization contexts include
53b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
54b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
55c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
56c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
57c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
589b94acceSBarry Smith 
593e081fefSLois Curfman McInnes    The user can open an alternative visualization context with
60b0a32e0cSBarry Smith    PetscViewerASCIIOpen() - output to a specified file.
619b94acceSBarry Smith 
6236851e7fSLois Curfman McInnes    Level: beginner
6336851e7fSLois Curfman McInnes 
649b94acceSBarry Smith .keywords: SNES, view
659b94acceSBarry Smith 
66b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
679b94acceSBarry Smith @*/
6863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer)
699b94acceSBarry Smith {
70fa9f3622SBarry Smith   SNESKSPEW           *kctx;
71dfbe8321SBarry Smith   PetscErrorCode      ierr;
7294b7f48cSBarry Smith   KSP                 ksp;
73a313700dSBarry Smith   const SNESType      type;
7432077d6dSBarry Smith   PetscTruth          iascii,isstring;
759b94acceSBarry Smith 
763a40ed3dSBarry Smith   PetscFunctionBegin;
774482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
783050cee2SBarry Smith   if (!viewer) {
797adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
803050cee2SBarry Smith   }
814482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
82c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,viewer,2);
8374679c65SBarry Smith 
8432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
85b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
8632077d6dSBarry Smith   if (iascii) {
877adad957SLisandro Dalcin     if (((PetscObject)snes)->prefix) {
887adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr);
893a7fca6bSBarry Smith     } else {
90b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
913a7fca6bSBarry Smith     }
92454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
93454a90a3SBarry Smith     if (type) {
94b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
95184914b5SBarry Smith     } else {
96b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
97184914b5SBarry Smith     }
98e7788613SBarry Smith     if (snes->ops->view) {
99b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
100e7788613SBarry Smith       ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);
101b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1020ef38995SBarry Smith     }
10377431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
104a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
10570441072SBarry Smith                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
10677431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
10777431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
1089b94acceSBarry Smith     if (snes->ksp_ewconv) {
109fa9f3622SBarry Smith       kctx = (SNESKSPEW *)snes->kspconvctx;
1109b94acceSBarry Smith       if (kctx) {
11177431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
112a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
113a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
1149b94acceSBarry Smith       }
1159b94acceSBarry Smith     }
1160f5bd95cSBarry Smith   } else if (isstring) {
117454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
118b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
11919bcc07fSBarry Smith   }
12094b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
121b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
12294b7f48cSBarry Smith   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
123b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1243a40ed3dSBarry Smith   PetscFunctionReturn(0);
1259b94acceSBarry Smith }
1269b94acceSBarry Smith 
12776b2cf59SMatthew Knepley /*
12876b2cf59SMatthew Knepley   We retain a list of functions that also take SNES command
12976b2cf59SMatthew Knepley   line options. These are called at the end SNESSetFromOptions()
13076b2cf59SMatthew Knepley */
13176b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5
132a7cc72afSBarry Smith static PetscInt numberofsetfromoptions;
1336849ba73SBarry Smith static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
13476b2cf59SMatthew Knepley 
135e74ef692SMatthew Knepley #undef __FUNCT__
136e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker"
137ac226902SBarry Smith /*@C
13876b2cf59SMatthew Knepley   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
13976b2cf59SMatthew Knepley 
14076b2cf59SMatthew Knepley   Not Collective
14176b2cf59SMatthew Knepley 
14276b2cf59SMatthew Knepley   Input Parameter:
14376b2cf59SMatthew Knepley . snescheck - function that checks for options
14476b2cf59SMatthew Knepley 
14576b2cf59SMatthew Knepley   Level: developer
14676b2cf59SMatthew Knepley 
14776b2cf59SMatthew Knepley .seealso: SNESSetFromOptions()
14876b2cf59SMatthew Knepley @*/
14963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
15076b2cf59SMatthew Knepley {
15176b2cf59SMatthew Knepley   PetscFunctionBegin;
15276b2cf59SMatthew Knepley   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
15377431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
15476b2cf59SMatthew Knepley   }
15576b2cf59SMatthew Knepley   othersetfromoptions[numberofsetfromoptions++] = snescheck;
15676b2cf59SMatthew Knepley   PetscFunctionReturn(0);
15776b2cf59SMatthew Knepley }
15876b2cf59SMatthew Knepley 
1594a2ae208SSatish Balay #undef __FUNCT__
1604a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions"
1619b94acceSBarry Smith /*@
16294b7f48cSBarry Smith    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
1639b94acceSBarry Smith 
164c7afd0dbSLois Curfman McInnes    Collective on SNES
165c7afd0dbSLois Curfman McInnes 
1669b94acceSBarry Smith    Input Parameter:
1679b94acceSBarry Smith .  snes - the SNES context
1689b94acceSBarry Smith 
16936851e7fSLois Curfman McInnes    Options Database Keys:
1706831982aSBarry Smith +  -snes_type <type> - ls, tr, umls, umtr, test
17182738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
17282738288SBarry Smith                 of the change in the solution between steps
17370441072SBarry Smith .  -snes_atol <abstol> - absolute tolerance of residual norm
174b39c3a46SLois Curfman McInnes .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
175b39c3a46SLois Curfman McInnes .  -snes_max_it <max_it> - maximum number of iterations
176b39c3a46SLois Curfman McInnes .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
17750ffb88aSMatthew Knepley .  -snes_max_fail <max_fail> - maximum number of failures
178ddf469c8SBarry Smith .  -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops
179*a8054027SBarry Smith .  -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild)
180b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
1812492ecdbSBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear
18282738288SBarry Smith                                solver; hence iterations will continue until max_it
1831fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
18482738288SBarry Smith                                of convergence test
185e8105e01SRichard Katz .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
186e8105e01SRichard Katz                                        filename given prints to stdout
187a6570f20SBarry Smith .  -snes_monitor_solution - plots solution at each iteration
188a6570f20SBarry Smith .  -snes_monitor_residual - plots residual (not its norm) at each iteration
189a6570f20SBarry Smith .  -snes_monitor_solution_update - plots update to solution at each iteration
190a6570f20SBarry Smith .  -snes_monitor_draw - plots residual norm at each iteration
191e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
1925968eb51SBarry Smith .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
193fee2055bSBarry Smith -  -snes_converged_reason - print the reason for convergence/divergence after each solve
19482738288SBarry Smith 
19582738288SBarry Smith     Options Database for Eisenstat-Walker method:
196fa9f3622SBarry Smith +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
1974b27c08aSLois Curfman McInnes .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
19836851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
19936851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
20036851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
20136851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
20236851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
20336851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
20482738288SBarry Smith 
20511ca99fdSLois Curfman McInnes    Notes:
20611ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
20711ca99fdSLois Curfman McInnes    the users manual.
20883e2fdc7SBarry Smith 
20936851e7fSLois Curfman McInnes    Level: beginner
21036851e7fSLois Curfman McInnes 
2119b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
2129b94acceSBarry Smith 
21369ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
2149b94acceSBarry Smith @*/
21563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes)
2169b94acceSBarry Smith {
217f1af5d2fSBarry Smith   PetscTruth              flg;
218*a8054027SBarry Smith   PetscInt                i,indx,lag;
21985385478SLisandro Dalcin   const char              *deft = SNESLS;
22085385478SLisandro Dalcin   const char              *convtests[] = {"default","skip"};
22185385478SLisandro Dalcin   SNESKSPEW               *kctx = NULL;
222e8105e01SRichard Katz   char                    type[256], monfilename[PETSC_MAX_PATH_LEN];
22323d894e5SBarry Smith   PetscViewerASCIIMonitor monviewer;
22485385478SLisandro Dalcin   PetscErrorCode          ierr;
2259b94acceSBarry Smith 
2263a40ed3dSBarry Smith   PetscFunctionBegin;
2274482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
228ca161407SBarry Smith 
2297adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
230186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2317adad957SLisandro Dalcin     if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; }
232b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
233d64ed03dSBarry Smith     if (flg) {
234186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
2357adad957SLisandro Dalcin     } else if (!((PetscObject)snes)->type_name) {
236186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
237d64ed03dSBarry Smith     }
238909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
23993c39befSBarry Smith 
24087828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
24170441072SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
242186905e3SBarry Smith 
24387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
244b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
245b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
24650ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
247ddf469c8SBarry Smith     ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr);
24885385478SLisandro Dalcin 
249*a8054027SBarry Smith     ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr);
250*a8054027SBarry Smith     if (flg) {
251*a8054027SBarry Smith       ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr);
252*a8054027SBarry Smith     }
253*a8054027SBarry Smith 
25485385478SLisandro Dalcin     ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr);
25585385478SLisandro Dalcin     if (flg) {
25685385478SLisandro Dalcin       switch (indx) {
2577f7931b9SBarry Smith       case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break;
2587f7931b9SBarry Smith       case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);    break;
25985385478SLisandro Dalcin       }
26085385478SLisandro Dalcin     }
26185385478SLisandro Dalcin 
2625968eb51SBarry Smith     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
2635968eb51SBarry Smith     if (flg) {
2645968eb51SBarry Smith       snes->printreason = PETSC_TRUE;
2655968eb51SBarry Smith     }
266186905e3SBarry Smith 
26785385478SLisandro Dalcin     kctx = (SNESKSPEW *)snes->kspconvctx;
26885385478SLisandro Dalcin 
269fa9f3622SBarry Smith     ierr = PetscOptionsTruth("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr);
270186905e3SBarry Smith 
271fa9f3622SBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
272fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
273fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
274fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
275fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
276fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
277fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
278186905e3SBarry Smith 
279a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",&flg);CHKERRQ(ierr);
280a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);}
281eabae89aSBarry Smith 
282a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
283e8105e01SRichard Katz     if (flg) {
2847adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
28523d894e5SBarry Smith       ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
286e8105e01SRichard Katz     }
287eabae89aSBarry Smith 
288a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
289eabae89aSBarry Smith     if (flg) {
2907adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
291f1bef1bcSMatthew Knepley       ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr);
292e8105e01SRichard Katz     }
293eabae89aSBarry Smith 
294a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
295eabae89aSBarry Smith     if (flg) {
2967adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
29723d894e5SBarry Smith       ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
298eabae89aSBarry Smith     }
299eabae89aSBarry Smith 
300a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",&flg);CHKERRQ(ierr);
301a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);}
302a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",&flg);CHKERRQ(ierr);
303a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);}
304a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",&flg);CHKERRQ(ierr);
305a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);}
306a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",&flg);CHKERRQ(ierr);
307a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
308e24b481bSBarry Smith 
309b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
3104b27c08aSLois Curfman McInnes     if (flg) {
311186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
312ae15b995SBarry Smith       ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr);
3139b94acceSBarry Smith     }
314639f9d9dSBarry Smith 
31576b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
31676b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
31776b2cf59SMatthew Knepley     }
31876b2cf59SMatthew Knepley 
319e7788613SBarry Smith     if (snes->ops->setfromoptions) {
320e7788613SBarry Smith       ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr);
321639f9d9dSBarry Smith     }
322b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
3234bbc92c1SBarry Smith 
32485385478SLisandro Dalcin   ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr);
32593993e2dSLois Curfman McInnes 
3263a40ed3dSBarry Smith   PetscFunctionReturn(0);
3279b94acceSBarry Smith }
3289b94acceSBarry Smith 
329a847f771SSatish Balay 
3304a2ae208SSatish Balay #undef __FUNCT__
3314a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
3329b94acceSBarry Smith /*@
3339b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3349b94acceSBarry Smith    the nonlinear solvers.
3359b94acceSBarry Smith 
336fee21e36SBarry Smith    Collective on SNES
337fee21e36SBarry Smith 
338c7afd0dbSLois Curfman McInnes    Input Parameters:
339c7afd0dbSLois Curfman McInnes +  snes - the SNES context
340c7afd0dbSLois Curfman McInnes -  usrP - optional user context
341c7afd0dbSLois Curfman McInnes 
34236851e7fSLois Curfman McInnes    Level: intermediate
34336851e7fSLois Curfman McInnes 
3449b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3459b94acceSBarry Smith 
3469b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3479b94acceSBarry Smith @*/
34863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP)
3499b94acceSBarry Smith {
3503a40ed3dSBarry Smith   PetscFunctionBegin;
3514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3529b94acceSBarry Smith   snes->user		= usrP;
3533a40ed3dSBarry Smith   PetscFunctionReturn(0);
3549b94acceSBarry Smith }
35574679c65SBarry Smith 
3564a2ae208SSatish Balay #undef __FUNCT__
3574a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3589b94acceSBarry Smith /*@C
3599b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3609b94acceSBarry Smith    nonlinear solvers.
3619b94acceSBarry Smith 
362c7afd0dbSLois Curfman McInnes    Not Collective
363c7afd0dbSLois Curfman McInnes 
3649b94acceSBarry Smith    Input Parameter:
3659b94acceSBarry Smith .  snes - SNES context
3669b94acceSBarry Smith 
3679b94acceSBarry Smith    Output Parameter:
3689b94acceSBarry Smith .  usrP - user context
3699b94acceSBarry Smith 
37036851e7fSLois Curfman McInnes    Level: intermediate
37136851e7fSLois Curfman McInnes 
3729b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3739b94acceSBarry Smith 
3749b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3759b94acceSBarry Smith @*/
37663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP)
3779b94acceSBarry Smith {
3783a40ed3dSBarry Smith   PetscFunctionBegin;
3794482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3809b94acceSBarry Smith   *usrP = snes->user;
3813a40ed3dSBarry Smith   PetscFunctionReturn(0);
3829b94acceSBarry Smith }
38374679c65SBarry Smith 
3844a2ae208SSatish Balay #undef __FUNCT__
3854a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3869b94acceSBarry Smith /*@
387c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
388c8228a4eSBarry Smith    at this time.
3899b94acceSBarry Smith 
390c7afd0dbSLois Curfman McInnes    Not Collective
391c7afd0dbSLois Curfman McInnes 
3929b94acceSBarry Smith    Input Parameter:
3939b94acceSBarry Smith .  snes - SNES context
3949b94acceSBarry Smith 
3959b94acceSBarry Smith    Output Parameter:
3969b94acceSBarry Smith .  iter - iteration number
3979b94acceSBarry Smith 
398c8228a4eSBarry Smith    Notes:
399c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
400c8228a4eSBarry Smith 
401c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
40208405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
40308405cd6SLois Curfman McInnes .vb
40408405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
40508405cd6SLois Curfman McInnes       if (!(it % 2)) {
40608405cd6SLois Curfman McInnes         [compute Jacobian here]
40708405cd6SLois Curfman McInnes       }
40808405cd6SLois Curfman McInnes .ve
409c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
41008405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
411c8228a4eSBarry Smith 
41236851e7fSLois Curfman McInnes    Level: intermediate
41336851e7fSLois Curfman McInnes 
4142b668275SBarry Smith .keywords: SNES, nonlinear, get, iteration, number,
4152b668275SBarry Smith 
416b850b91aSLisandro Dalcin .seealso:   SNESGetFunctionNorm(), SNESGetLinearSolveIterations()
4179b94acceSBarry Smith @*/
41863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter)
4199b94acceSBarry Smith {
4203a40ed3dSBarry Smith   PetscFunctionBegin;
4214482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4224482741eSBarry Smith   PetscValidIntPointer(iter,2);
4239b94acceSBarry Smith   *iter = snes->iter;
4243a40ed3dSBarry Smith   PetscFunctionReturn(0);
4259b94acceSBarry Smith }
42674679c65SBarry Smith 
4274a2ae208SSatish Balay #undef __FUNCT__
4284a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
4299b94acceSBarry Smith /*@
4309b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
4319b94acceSBarry Smith    with SNESSSetFunction().
4329b94acceSBarry Smith 
433c7afd0dbSLois Curfman McInnes    Collective on SNES
434c7afd0dbSLois Curfman McInnes 
4359b94acceSBarry Smith    Input Parameter:
4369b94acceSBarry Smith .  snes - SNES context
4379b94acceSBarry Smith 
4389b94acceSBarry Smith    Output Parameter:
4399b94acceSBarry Smith .  fnorm - 2-norm of function
4409b94acceSBarry Smith 
44136851e7fSLois Curfman McInnes    Level: intermediate
44236851e7fSLois Curfman McInnes 
4439b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
444a86d99e1SLois Curfman McInnes 
445b850b91aSLisandro Dalcin .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations()
4469b94acceSBarry Smith @*/
44771f87433Sdalcinl PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscReal *fnorm)
4489b94acceSBarry Smith {
4493a40ed3dSBarry Smith   PetscFunctionBegin;
4504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4514482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
4529b94acceSBarry Smith   *fnorm = snes->norm;
4533a40ed3dSBarry Smith   PetscFunctionReturn(0);
4549b94acceSBarry Smith }
45574679c65SBarry Smith 
4564a2ae208SSatish Balay #undef __FUNCT__
457b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetNonlinearStepFailures"
4589b94acceSBarry Smith /*@
459b850b91aSLisandro Dalcin    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
4609b94acceSBarry Smith    attempted by the nonlinear solver.
4619b94acceSBarry Smith 
462c7afd0dbSLois Curfman McInnes    Not Collective
463c7afd0dbSLois Curfman McInnes 
4649b94acceSBarry Smith    Input Parameter:
4659b94acceSBarry Smith .  snes - SNES context
4669b94acceSBarry Smith 
4679b94acceSBarry Smith    Output Parameter:
4689b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4699b94acceSBarry Smith 
470c96a6f78SLois Curfman McInnes    Notes:
471c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
472c96a6f78SLois Curfman McInnes 
47336851e7fSLois Curfman McInnes    Level: intermediate
47436851e7fSLois Curfman McInnes 
4759b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
47658ebbce7SBarry Smith 
477e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
47858ebbce7SBarry Smith           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
4799b94acceSBarry Smith @*/
480b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails)
4819b94acceSBarry Smith {
4823a40ed3dSBarry Smith   PetscFunctionBegin;
4834482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4844482741eSBarry Smith   PetscValidIntPointer(nfails,2);
48550ffb88aSMatthew Knepley   *nfails = snes->numFailures;
48650ffb88aSMatthew Knepley   PetscFunctionReturn(0);
48750ffb88aSMatthew Knepley }
48850ffb88aSMatthew Knepley 
48950ffb88aSMatthew Knepley #undef __FUNCT__
490b850b91aSLisandro Dalcin #define __FUNCT__ "SNESSetMaxNonlinearStepFailures"
49150ffb88aSMatthew Knepley /*@
492b850b91aSLisandro Dalcin    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
49350ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
49450ffb88aSMatthew Knepley 
49550ffb88aSMatthew Knepley    Not Collective
49650ffb88aSMatthew Knepley 
49750ffb88aSMatthew Knepley    Input Parameters:
49850ffb88aSMatthew Knepley +  snes     - SNES context
49950ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
50050ffb88aSMatthew Knepley 
50150ffb88aSMatthew Knepley    Level: intermediate
50250ffb88aSMatthew Knepley 
50350ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
50458ebbce7SBarry Smith 
505e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
50658ebbce7SBarry Smith           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
50750ffb88aSMatthew Knepley @*/
508b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
50950ffb88aSMatthew Knepley {
51050ffb88aSMatthew Knepley   PetscFunctionBegin;
5114482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
51250ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
51350ffb88aSMatthew Knepley   PetscFunctionReturn(0);
51450ffb88aSMatthew Knepley }
51550ffb88aSMatthew Knepley 
51650ffb88aSMatthew Knepley #undef __FUNCT__
517b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetMaxNonlinearStepFailures"
51850ffb88aSMatthew Knepley /*@
519b850b91aSLisandro Dalcin    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
52050ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
52150ffb88aSMatthew Knepley 
52250ffb88aSMatthew Knepley    Not Collective
52350ffb88aSMatthew Knepley 
52450ffb88aSMatthew Knepley    Input Parameter:
52550ffb88aSMatthew Knepley .  snes     - SNES context
52650ffb88aSMatthew Knepley 
52750ffb88aSMatthew Knepley    Output Parameter:
52850ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
52950ffb88aSMatthew Knepley 
53050ffb88aSMatthew Knepley    Level: intermediate
53150ffb88aSMatthew Knepley 
53250ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
53358ebbce7SBarry Smith 
534e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
53558ebbce7SBarry Smith           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
53658ebbce7SBarry Smith 
53750ffb88aSMatthew Knepley @*/
538b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
53950ffb88aSMatthew Knepley {
54050ffb88aSMatthew Knepley   PetscFunctionBegin;
5414482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5424482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
54350ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
5443a40ed3dSBarry Smith   PetscFunctionReturn(0);
5459b94acceSBarry Smith }
546a847f771SSatish Balay 
5474a2ae208SSatish Balay #undef __FUNCT__
5482541af92SBarry Smith #define __FUNCT__ "SNESGetNumberFunctionEvals"
5492541af92SBarry Smith /*@
5502541af92SBarry Smith    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
5512541af92SBarry Smith      done by SNES.
5522541af92SBarry Smith 
5532541af92SBarry Smith    Not Collective
5542541af92SBarry Smith 
5552541af92SBarry Smith    Input Parameter:
5562541af92SBarry Smith .  snes     - SNES context
5572541af92SBarry Smith 
5582541af92SBarry Smith    Output Parameter:
5592541af92SBarry Smith .  nfuncs - number of evaluations
5602541af92SBarry Smith 
5612541af92SBarry Smith    Level: intermediate
5622541af92SBarry Smith 
5632541af92SBarry Smith .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
56458ebbce7SBarry Smith 
565e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures()
5662541af92SBarry Smith @*/
5672541af92SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
5682541af92SBarry Smith {
5692541af92SBarry Smith   PetscFunctionBegin;
5702541af92SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5712541af92SBarry Smith   PetscValidIntPointer(nfuncs,2);
5722541af92SBarry Smith   *nfuncs = snes->nfuncs;
5732541af92SBarry Smith   PetscFunctionReturn(0);
5742541af92SBarry Smith }
5752541af92SBarry Smith 
5762541af92SBarry Smith #undef __FUNCT__
5773d4c4710SBarry Smith #define __FUNCT__ "SNESGetLinearSolveFailures"
5783d4c4710SBarry Smith /*@
5793d4c4710SBarry Smith    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
5803d4c4710SBarry Smith    linear solvers.
5813d4c4710SBarry Smith 
5823d4c4710SBarry Smith    Not Collective
5833d4c4710SBarry Smith 
5843d4c4710SBarry Smith    Input Parameter:
5853d4c4710SBarry Smith .  snes - SNES context
5863d4c4710SBarry Smith 
5873d4c4710SBarry Smith    Output Parameter:
5883d4c4710SBarry Smith .  nfails - number of failed solves
5893d4c4710SBarry Smith 
5903d4c4710SBarry Smith    Notes:
5913d4c4710SBarry Smith    This counter is reset to zero for each successive call to SNESSolve().
5923d4c4710SBarry Smith 
5933d4c4710SBarry Smith    Level: intermediate
5943d4c4710SBarry Smith 
5953d4c4710SBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
59658ebbce7SBarry Smith 
597e1c61ce8SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures()
5983d4c4710SBarry Smith @*/
5993d4c4710SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails)
6003d4c4710SBarry Smith {
6013d4c4710SBarry Smith   PetscFunctionBegin;
6023d4c4710SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6033d4c4710SBarry Smith   PetscValidIntPointer(nfails,2);
6043d4c4710SBarry Smith   *nfails = snes->numLinearSolveFailures;
6053d4c4710SBarry Smith   PetscFunctionReturn(0);
6063d4c4710SBarry Smith }
6073d4c4710SBarry Smith 
6083d4c4710SBarry Smith #undef __FUNCT__
6093d4c4710SBarry Smith #define __FUNCT__ "SNESSetMaxLinearSolveFailures"
6103d4c4710SBarry Smith /*@
6113d4c4710SBarry Smith    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
6123d4c4710SBarry Smith    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
6133d4c4710SBarry Smith 
6143d4c4710SBarry Smith    Collective on SNES
6153d4c4710SBarry Smith 
6163d4c4710SBarry Smith    Input Parameters:
6173d4c4710SBarry Smith +  snes     - SNES context
6183d4c4710SBarry Smith -  maxFails - maximum allowed linear solve failures
6193d4c4710SBarry Smith 
6203d4c4710SBarry Smith    Level: intermediate
6213d4c4710SBarry Smith 
622a6796414SBarry Smith    Notes: By default this is 0; that is SNES returns on the first failed linear solve
6233d4c4710SBarry Smith 
6243d4c4710SBarry Smith .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
6253d4c4710SBarry Smith 
62658ebbce7SBarry Smith .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
6273d4c4710SBarry Smith @*/
6283d4c4710SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
6293d4c4710SBarry Smith {
6303d4c4710SBarry Smith   PetscFunctionBegin;
6313d4c4710SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6323d4c4710SBarry Smith   snes->maxLinearSolveFailures = maxFails;
6333d4c4710SBarry Smith   PetscFunctionReturn(0);
6343d4c4710SBarry Smith }
6353d4c4710SBarry Smith 
6363d4c4710SBarry Smith #undef __FUNCT__
6373d4c4710SBarry Smith #define __FUNCT__ "SNESGetMaxLinearSolveFailures"
6383d4c4710SBarry Smith /*@
6393d4c4710SBarry Smith    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
6403d4c4710SBarry Smith      are allowed before SNES terminates
6413d4c4710SBarry Smith 
6423d4c4710SBarry Smith    Not Collective
6433d4c4710SBarry Smith 
6443d4c4710SBarry Smith    Input Parameter:
6453d4c4710SBarry Smith .  snes     - SNES context
6463d4c4710SBarry Smith 
6473d4c4710SBarry Smith    Output Parameter:
6483d4c4710SBarry Smith .  maxFails - maximum of unsuccessful solves allowed
6493d4c4710SBarry Smith 
6503d4c4710SBarry Smith    Level: intermediate
6513d4c4710SBarry Smith 
6523d4c4710SBarry Smith    Notes: By default this is 1; that is SNES returns on the first failed linear solve
6533d4c4710SBarry Smith 
6543d4c4710SBarry Smith .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
6553d4c4710SBarry Smith 
656e1c61ce8SBarry Smith .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(),
6573d4c4710SBarry Smith @*/
6583d4c4710SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
6593d4c4710SBarry Smith {
6603d4c4710SBarry Smith   PetscFunctionBegin;
6613d4c4710SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6623d4c4710SBarry Smith   PetscValidIntPointer(maxFails,2);
6633d4c4710SBarry Smith   *maxFails = snes->maxLinearSolveFailures;
6643d4c4710SBarry Smith   PetscFunctionReturn(0);
6653d4c4710SBarry Smith }
6663d4c4710SBarry Smith 
6673d4c4710SBarry Smith #undef __FUNCT__
668b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetLinearSolveIterations"
669c96a6f78SLois Curfman McInnes /*@
670b850b91aSLisandro Dalcin    SNESGetLinearSolveIterations - Gets the total number of linear iterations
671c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
672c96a6f78SLois Curfman McInnes 
673c7afd0dbSLois Curfman McInnes    Not Collective
674c7afd0dbSLois Curfman McInnes 
675c96a6f78SLois Curfman McInnes    Input Parameter:
676c96a6f78SLois Curfman McInnes .  snes - SNES context
677c96a6f78SLois Curfman McInnes 
678c96a6f78SLois Curfman McInnes    Output Parameter:
679c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
680c96a6f78SLois Curfman McInnes 
681c96a6f78SLois Curfman McInnes    Notes:
682c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
683c96a6f78SLois Curfman McInnes 
68436851e7fSLois Curfman McInnes    Level: intermediate
68536851e7fSLois Curfman McInnes 
686c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
6872b668275SBarry Smith 
68858ebbce7SBarry Smith .seealso:  SNESGetIterationNumber(), SNESGetFunctionNorm()S, NESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures()
689c96a6f78SLois Curfman McInnes @*/
690b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveIterations(SNES snes,PetscInt* lits)
691c96a6f78SLois Curfman McInnes {
6923a40ed3dSBarry Smith   PetscFunctionBegin;
6934482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6944482741eSBarry Smith   PetscValidIntPointer(lits,2);
695c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
6963a40ed3dSBarry Smith   PetscFunctionReturn(0);
697c96a6f78SLois Curfman McInnes }
698c96a6f78SLois Curfman McInnes 
6994a2ae208SSatish Balay #undef __FUNCT__
70094b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
70152baeb72SSatish Balay /*@
70294b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
7039b94acceSBarry Smith 
70494b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
705c7afd0dbSLois Curfman McInnes 
7069b94acceSBarry Smith    Input Parameter:
7079b94acceSBarry Smith .  snes - the SNES context
7089b94acceSBarry Smith 
7099b94acceSBarry Smith    Output Parameter:
71094b7f48cSBarry Smith .  ksp - the KSP context
7119b94acceSBarry Smith 
7129b94acceSBarry Smith    Notes:
71394b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
7149b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
7152999313aSBarry Smith    PC contexts as well.
7169b94acceSBarry Smith 
71736851e7fSLois Curfman McInnes    Level: beginner
71836851e7fSLois Curfman McInnes 
71994b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
7209b94acceSBarry Smith 
7212999313aSBarry Smith .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
7229b94acceSBarry Smith @*/
72363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp)
7249b94acceSBarry Smith {
7253a40ed3dSBarry Smith   PetscFunctionBegin;
7264482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7274482741eSBarry Smith   PetscValidPointer(ksp,2);
72894b7f48cSBarry Smith   *ksp = snes->ksp;
7293a40ed3dSBarry Smith   PetscFunctionReturn(0);
7309b94acceSBarry Smith }
73182bf6240SBarry Smith 
7324a2ae208SSatish Balay #undef __FUNCT__
7332999313aSBarry Smith #define __FUNCT__ "SNESSetKSP"
7342999313aSBarry Smith /*@
7352999313aSBarry Smith    SNESSetKSP - Sets a KSP context for the SNES object to use
7362999313aSBarry Smith 
7372999313aSBarry Smith    Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
7382999313aSBarry Smith 
7392999313aSBarry Smith    Input Parameters:
7402999313aSBarry Smith +  snes - the SNES context
7412999313aSBarry Smith -  ksp - the KSP context
7422999313aSBarry Smith 
7432999313aSBarry Smith    Notes:
7442999313aSBarry Smith    The SNES object already has its KSP object, you can obtain with SNESGetKSP()
7452999313aSBarry Smith    so this routine is rarely needed.
7462999313aSBarry Smith 
7472999313aSBarry Smith    The KSP object that is already in the SNES object has its reference count
7482999313aSBarry Smith    decreased by one.
7492999313aSBarry Smith 
7502999313aSBarry Smith    Level: developer
7512999313aSBarry Smith 
7522999313aSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
7532999313aSBarry Smith 
7542999313aSBarry Smith .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
7552999313aSBarry Smith @*/
7562999313aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetKSP(SNES snes,KSP ksp)
7572999313aSBarry Smith {
7582999313aSBarry Smith   PetscErrorCode ierr;
7592999313aSBarry Smith 
7602999313aSBarry Smith   PetscFunctionBegin;
7612999313aSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7622999313aSBarry Smith   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
7632999313aSBarry Smith   PetscCheckSameComm(snes,1,ksp,2);
7647dcf0eaaSdalcinl   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
765906ed7ccSBarry Smith   if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);}
7662999313aSBarry Smith   snes->ksp = ksp;
7672999313aSBarry Smith   PetscFunctionReturn(0);
7682999313aSBarry Smith }
7692999313aSBarry Smith 
7707adad957SLisandro Dalcin #if 0
7712999313aSBarry Smith #undef __FUNCT__
7724a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
7736849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
774e24b481bSBarry Smith {
775e24b481bSBarry Smith   PetscFunctionBegin;
776e24b481bSBarry Smith   PetscFunctionReturn(0);
777e24b481bSBarry Smith }
7787adad957SLisandro Dalcin #endif
779e24b481bSBarry Smith 
7809b94acceSBarry Smith /* -----------------------------------------------------------*/
7814a2ae208SSatish Balay #undef __FUNCT__
7824a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
78352baeb72SSatish Balay /*@
7849b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
7859b94acceSBarry Smith 
786c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
787c7afd0dbSLois Curfman McInnes 
788c7afd0dbSLois Curfman McInnes    Input Parameters:
789906ed7ccSBarry Smith .  comm - MPI communicator
7909b94acceSBarry Smith 
7919b94acceSBarry Smith    Output Parameter:
7929b94acceSBarry Smith .  outsnes - the new SNES context
7939b94acceSBarry Smith 
794c7afd0dbSLois Curfman McInnes    Options Database Keys:
795c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
796c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
797c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
798c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
799c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
800c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
801c1f60f51SBarry Smith 
80236851e7fSLois Curfman McInnes    Level: beginner
80336851e7fSLois Curfman McInnes 
8049b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
8059b94acceSBarry Smith 
806*a8054027SBarry Smith .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner()
807*a8054027SBarry Smith 
8089b94acceSBarry Smith @*/
80963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes)
8109b94acceSBarry Smith {
811dfbe8321SBarry Smith   PetscErrorCode      ierr;
8129b94acceSBarry Smith   SNES                snes;
813fa9f3622SBarry Smith   SNESKSPEW           *kctx;
81437fcc0dbSBarry Smith 
8153a40ed3dSBarry Smith   PetscFunctionBegin;
816ed1caa07SMatthew Knepley   PetscValidPointer(outsnes,2);
8178ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
8188ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
8198ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
8208ba1e511SMatthew Knepley #endif
8218ba1e511SMatthew Knepley 
822e7788613SBarry Smith   ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
8237adad957SLisandro Dalcin 
82485385478SLisandro Dalcin   snes->ops->converged    = SNESDefaultConverged;
8259b94acceSBarry Smith   snes->max_its           = 50;
8269750a799SBarry Smith   snes->max_funcs	  = 10000;
8279b94acceSBarry Smith   snes->norm		  = 0.0;
828b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
829b4874afaSBarry Smith   snes->ttol              = 0.0;
83070441072SBarry Smith   snes->abstol		  = 1.e-50;
8319b94acceSBarry Smith   snes->xtol		  = 1.e-8;
8324b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
8339b94acceSBarry Smith   snes->nfuncs            = 0;
83450ffb88aSMatthew Knepley   snes->numFailures       = 0;
83550ffb88aSMatthew Knepley   snes->maxFailures       = 1;
8367a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
837*a8054027SBarry Smith   snes->lagpreconditioner = 1;
838639f9d9dSBarry Smith   snes->numbermonitors    = 0;
8399b94acceSBarry Smith   snes->data              = 0;
8404dc4c822SBarry Smith   snes->setupcalled       = PETSC_FALSE;
841186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
8426f24a144SLois Curfman McInnes   snes->vwork             = 0;
8436f24a144SLois Curfman McInnes   snes->nwork             = 0;
844758f92a0SBarry Smith   snes->conv_hist_len     = 0;
845758f92a0SBarry Smith   snes->conv_hist_max     = 0;
846758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
847758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
848758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
849184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
8509b94acceSBarry Smith 
8513d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
8523d4c4710SBarry Smith   snes->maxLinearSolveFailures = 1;
8533d4c4710SBarry Smith 
8549b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
85538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr);
8569b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
8579b94acceSBarry Smith   kctx->version     = 2;
8589b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
8599b94acceSBarry Smith                              this was too large for some test cases */
8609b94acceSBarry Smith   kctx->rtol_last   = 0;
8619b94acceSBarry Smith   kctx->rtol_max    = .9;
8629b94acceSBarry Smith   kctx->gamma       = 1.0;
86371f87433Sdalcinl   kctx->alpha       = .5*(1.0 + sqrt(5.0));
86471f87433Sdalcinl   kctx->alpha2      = kctx->alpha;
8659b94acceSBarry Smith   kctx->threshold   = .1;
8669b94acceSBarry Smith   kctx->lresid_last = 0;
8679b94acceSBarry Smith   kctx->norm_last   = 0;
8689b94acceSBarry Smith 
86994b7f48cSBarry Smith   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
87052e6d16bSBarry Smith   ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
8715334005bSBarry Smith 
8729b94acceSBarry Smith   *outsnes = snes;
87300036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
8743a40ed3dSBarry Smith   PetscFunctionReturn(0);
8759b94acceSBarry Smith }
8769b94acceSBarry Smith 
8774a2ae208SSatish Balay #undef __FUNCT__
8784a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
8799b94acceSBarry Smith /*@C
8809b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
8819b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
8829b94acceSBarry Smith    equations.
8839b94acceSBarry Smith 
884fee21e36SBarry Smith    Collective on SNES
885fee21e36SBarry Smith 
886c7afd0dbSLois Curfman McInnes    Input Parameters:
887c7afd0dbSLois Curfman McInnes +  snes - the SNES context
888c7afd0dbSLois Curfman McInnes .  r - vector to store function value
889de044059SHong Zhang .  func - function evaluation routine
890c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
891c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
8929b94acceSBarry Smith 
893c7afd0dbSLois Curfman McInnes    Calling sequence of func:
8948d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
895c7afd0dbSLois Curfman McInnes 
896313e4042SLois Curfman McInnes .  f - function vector
897c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
8989b94acceSBarry Smith 
8999b94acceSBarry Smith    Notes:
9009b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
9019b94acceSBarry Smith $      f'(x) x = -f(x),
902c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
9039b94acceSBarry Smith 
90436851e7fSLois Curfman McInnes    Level: beginner
90536851e7fSLois Curfman McInnes 
9069b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
9079b94acceSBarry Smith 
908a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
9099b94acceSBarry Smith @*/
91063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
9119b94acceSBarry Smith {
91285385478SLisandro Dalcin   PetscErrorCode ierr;
9133a40ed3dSBarry Smith   PetscFunctionBegin;
9144482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9154482741eSBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
916c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,r,2);
91785385478SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr);
91885385478SLisandro Dalcin   if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); }
919e7788613SBarry Smith   snes->ops->computefunction = func;
92085385478SLisandro Dalcin   snes->vec_func             = r;
9219b94acceSBarry Smith   snes->funP                 = ctx;
9223a40ed3dSBarry Smith   PetscFunctionReturn(0);
9239b94acceSBarry Smith }
9249b94acceSBarry Smith 
9253ab0aad5SBarry Smith /* --------------------------------------------------------------- */
9263ab0aad5SBarry Smith #undef __FUNCT__
9271096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs"
9281096aae1SMatthew Knepley /*@C
9291096aae1SMatthew Knepley    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
9301096aae1SMatthew Knepley    it assumes a zero right hand side.
9311096aae1SMatthew Knepley 
9321096aae1SMatthew Knepley    Collective on SNES
9331096aae1SMatthew Knepley 
9341096aae1SMatthew Knepley    Input Parameter:
9351096aae1SMatthew Knepley .  snes - the SNES context
9361096aae1SMatthew Knepley 
9371096aae1SMatthew Knepley    Output Parameter:
938bc08b0f1SBarry Smith .  rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null
9391096aae1SMatthew Knepley 
9401096aae1SMatthew Knepley    Level: intermediate
9411096aae1SMatthew Knepley 
9421096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side
9431096aae1SMatthew Knepley 
94485385478SLisandro Dalcin .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
9451096aae1SMatthew Knepley @*/
9461096aae1SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs)
9471096aae1SMatthew Knepley {
9481096aae1SMatthew Knepley   PetscFunctionBegin;
9491096aae1SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9501096aae1SMatthew Knepley   PetscValidPointer(rhs,2);
95185385478SLisandro Dalcin   *rhs = snes->vec_rhs;
9521096aae1SMatthew Knepley   PetscFunctionReturn(0);
9531096aae1SMatthew Knepley }
9541096aae1SMatthew Knepley 
9551096aae1SMatthew Knepley #undef __FUNCT__
9564a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
9579b94acceSBarry Smith /*@
95836851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
9599b94acceSBarry Smith                          SNESSetFunction().
9609b94acceSBarry Smith 
961c7afd0dbSLois Curfman McInnes    Collective on SNES
962c7afd0dbSLois Curfman McInnes 
9639b94acceSBarry Smith    Input Parameters:
964c7afd0dbSLois Curfman McInnes +  snes - the SNES context
965c7afd0dbSLois Curfman McInnes -  x - input vector
9669b94acceSBarry Smith 
9679b94acceSBarry Smith    Output Parameter:
9683638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
9699b94acceSBarry Smith 
9701bffabb2SLois Curfman McInnes    Notes:
97136851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
97236851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
97336851e7fSLois Curfman McInnes    themselves.
97436851e7fSLois Curfman McInnes 
97536851e7fSLois Curfman McInnes    Level: developer
97636851e7fSLois Curfman McInnes 
9779b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
9789b94acceSBarry Smith 
979a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
9809b94acceSBarry Smith @*/
98163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y)
9829b94acceSBarry Smith {
983dfbe8321SBarry Smith   PetscErrorCode ierr;
9849b94acceSBarry Smith 
9853a40ed3dSBarry Smith   PetscFunctionBegin;
9864482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9874482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
9884482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
989c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
990c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
991184914b5SBarry Smith 
992d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
993e7788613SBarry Smith   if (snes->ops->computefunction) {
994d64ed03dSBarry Smith     PetscStackPush("SNES user function");
995e9a2bbcdSBarry Smith     CHKMEMQ;
996e7788613SBarry Smith     ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);
997e9a2bbcdSBarry Smith     CHKMEMQ;
998d64ed03dSBarry Smith     PetscStackPop;
999d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
100019717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
100119717074SBarry Smith     }
1002d5e45103SBarry Smith     CHKERRQ(ierr);
100385385478SLisandro Dalcin   } else if (snes->vec_rhs) {
10041096aae1SMatthew Knepley     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
10051096aae1SMatthew Knepley   } else {
10061096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
10071096aae1SMatthew Knepley   }
100885385478SLisandro Dalcin   if (snes->vec_rhs) {
100985385478SLisandro Dalcin     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
10103ab0aad5SBarry Smith   }
1011ae3c334cSLois Curfman McInnes   snes->nfuncs++;
1012d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
10133a40ed3dSBarry Smith   PetscFunctionReturn(0);
10149b94acceSBarry Smith }
10159b94acceSBarry Smith 
10164a2ae208SSatish Balay #undef __FUNCT__
10174a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
101862fef451SLois Curfman McInnes /*@
101962fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
102062fef451SLois Curfman McInnes    set with SNESSetJacobian().
102162fef451SLois Curfman McInnes 
1022c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1023c7afd0dbSLois Curfman McInnes 
102462fef451SLois Curfman McInnes    Input Parameters:
1025c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1026c7afd0dbSLois Curfman McInnes -  x - input vector
102762fef451SLois Curfman McInnes 
102862fef451SLois Curfman McInnes    Output Parameters:
1029c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
103062fef451SLois Curfman McInnes .  B - optional preconditioning matrix
10312b668275SBarry Smith -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
1032fee21e36SBarry Smith 
103362fef451SLois Curfman McInnes    Notes:
103462fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
103562fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
103662fef451SLois Curfman McInnes 
103794b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
1038dc5a77f8SLois Curfman McInnes    flag parameter.
103962fef451SLois Curfman McInnes 
104036851e7fSLois Curfman McInnes    Level: developer
104136851e7fSLois Curfman McInnes 
104262fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
104362fef451SLois Curfman McInnes 
10442b668275SBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure
104562fef451SLois Curfman McInnes @*/
104663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
10479b94acceSBarry Smith {
1048dfbe8321SBarry Smith   PetscErrorCode ierr;
10493a40ed3dSBarry Smith 
10503a40ed3dSBarry Smith   PetscFunctionBegin;
10514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10524482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
10534482741eSBarry Smith   PetscValidPointer(flg,5);
1054c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
1055e7788613SBarry Smith   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1056d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1057c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
1058d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
1059dc67937bSBarry Smith   CHKMEMQ;
1060e7788613SBarry Smith   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1061dc67937bSBarry Smith   CHKMEMQ;
1062d64ed03dSBarry Smith   PetscStackPop;
1063d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1064*a8054027SBarry Smith 
1065*a8054027SBarry Smith   if (snes->lagpreconditioner == -1) {
1066*a8054027SBarry Smith     *flg = SAME_PRECONDITIONER;
1067*a8054027SBarry Smith     ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr);
1068*a8054027SBarry Smith   } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) {
1069*a8054027SBarry Smith     *flg = SAME_PRECONDITIONER;
1070*a8054027SBarry Smith     ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr);
1071*a8054027SBarry Smith   }
1072*a8054027SBarry Smith 
10736d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
10746ce558aeSBarry Smith   /* PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
10756ce558aeSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,4);   */
10763a40ed3dSBarry Smith   PetscFunctionReturn(0);
10779b94acceSBarry Smith }
10789b94acceSBarry Smith 
10794a2ae208SSatish Balay #undef __FUNCT__
10804a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
10819b94acceSBarry Smith /*@C
10829b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1083044dda88SLois Curfman McInnes    location to store the matrix.
10849b94acceSBarry Smith 
1085c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1086c7afd0dbSLois Curfman McInnes 
10879b94acceSBarry Smith    Input Parameters:
1088c7afd0dbSLois Curfman McInnes +  snes - the SNES context
10899b94acceSBarry Smith .  A - Jacobian matrix
10909b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
10919b94acceSBarry Smith .  func - Jacobian evaluation routine
1092c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
10932cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
10949b94acceSBarry Smith 
10959b94acceSBarry Smith    Calling sequence of func:
10968d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10979b94acceSBarry Smith 
1098c7afd0dbSLois Curfman McInnes +  x - input vector
10999b94acceSBarry Smith .  A - Jacobian matrix
11009b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1101ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
11022b668275SBarry Smith    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1103c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
11049b94acceSBarry Smith 
11059b94acceSBarry Smith    Notes:
110694b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
11072cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1108ac21db08SLois Curfman McInnes 
1109ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
11109b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
11119b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
11129b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
11139b94acceSBarry Smith    throughout the global iterations.
11149b94acceSBarry Smith 
111516913363SBarry Smith    If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on
111616913363SBarry Smith    each matrix.
111716913363SBarry Smith 
111836851e7fSLois Curfman McInnes    Level: beginner
111936851e7fSLois Curfman McInnes 
11209b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
11219b94acceSBarry Smith 
11223ec795f1SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
11239b94acceSBarry Smith @*/
112463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
11259b94acceSBarry Smith {
1126dfbe8321SBarry Smith   PetscErrorCode ierr;
11273a7fca6bSBarry Smith 
11283a40ed3dSBarry Smith   PetscFunctionBegin;
11294482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
11304482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
11314482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
1132c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
1133c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
1134e7788613SBarry Smith   if (func) snes->ops->computejacobian = func;
11353a7fca6bSBarry Smith   if (ctx)  snes->jacP                 = ctx;
11363a7fca6bSBarry Smith   if (A) {
11377dcf0eaaSdalcinl     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
11383a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
11399b94acceSBarry Smith     snes->jacobian = A;
11403a7fca6bSBarry Smith   }
11413a7fca6bSBarry Smith   if (B) {
11427dcf0eaaSdalcinl     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
11433a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
11449b94acceSBarry Smith     snes->jacobian_pre = B;
11453a7fca6bSBarry Smith   }
11463a40ed3dSBarry Smith   PetscFunctionReturn(0);
11479b94acceSBarry Smith }
114862fef451SLois Curfman McInnes 
11494a2ae208SSatish Balay #undef __FUNCT__
11504a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
1151c2aafc4cSSatish Balay /*@C
1152b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1153b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1154b4fd4287SBarry Smith 
1155c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1156c7afd0dbSLois Curfman McInnes 
1157b4fd4287SBarry Smith    Input Parameter:
1158b4fd4287SBarry Smith .  snes - the nonlinear solver context
1159b4fd4287SBarry Smith 
1160b4fd4287SBarry Smith    Output Parameters:
1161c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1162b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
116370e92668SMatthew Knepley .  func - location to put Jacobian function (or PETSC_NULL)
116470e92668SMatthew Knepley -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1165fee21e36SBarry Smith 
116636851e7fSLois Curfman McInnes    Level: advanced
116736851e7fSLois Curfman McInnes 
1168b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1169b4fd4287SBarry Smith @*/
117070e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1171b4fd4287SBarry Smith {
11723a40ed3dSBarry Smith   PetscFunctionBegin;
11734482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1174b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
1175b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
1176e7788613SBarry Smith   if (func) *func = snes->ops->computejacobian;
117770e92668SMatthew Knepley   if (ctx)  *ctx  = snes->jacP;
11783a40ed3dSBarry Smith   PetscFunctionReturn(0);
1179b4fd4287SBarry Smith }
1180b4fd4287SBarry Smith 
11819b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
118263dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
11839b94acceSBarry Smith 
11844a2ae208SSatish Balay #undef __FUNCT__
11854a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
11869b94acceSBarry Smith /*@
11879b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1188272ac6f2SLois Curfman McInnes    of a nonlinear solver.
11899b94acceSBarry Smith 
1190fee21e36SBarry Smith    Collective on SNES
1191fee21e36SBarry Smith 
1192c7afd0dbSLois Curfman McInnes    Input Parameters:
119370e92668SMatthew Knepley .  snes - the SNES context
1194c7afd0dbSLois Curfman McInnes 
1195272ac6f2SLois Curfman McInnes    Notes:
1196272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1197272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1198272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1199272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1200272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1201272ac6f2SLois Curfman McInnes 
120236851e7fSLois Curfman McInnes    Level: advanced
120336851e7fSLois Curfman McInnes 
12049b94acceSBarry Smith .keywords: SNES, nonlinear, setup
12059b94acceSBarry Smith 
12069b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
12079b94acceSBarry Smith @*/
120870e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
12099b94acceSBarry Smith {
1210dfbe8321SBarry Smith   PetscErrorCode ierr;
121171f87433Sdalcinl   PetscTruth     flg;
12123a40ed3dSBarry Smith 
12133a40ed3dSBarry Smith   PetscFunctionBegin;
12144482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12154dc4c822SBarry Smith   if (snes->setupcalled) PetscFunctionReturn(0);
12169b94acceSBarry Smith 
12177adad957SLisandro Dalcin   if (!((PetscObject)snes)->type_name) {
121885385478SLisandro Dalcin     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
121985385478SLisandro Dalcin   }
122085385478SLisandro Dalcin 
12217adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1222c1f60f51SBarry Smith   /*
1223c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1224dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1225c1f60f51SBarry Smith   */
1226c1f60f51SBarry Smith   if (flg) {
1227c1f60f51SBarry Smith     Mat J;
1228fef1beadSBarry Smith     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
12293ec795f1SBarry Smith     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1230ae15b995SBarry Smith     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
12313a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
12323a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1233c1f60f51SBarry Smith   }
123445fc7adcSBarry Smith 
123503c60df9SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT)
12367adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
123745fc7adcSBarry Smith   if (flg) {
123845fc7adcSBarry Smith     Mat J;
123945fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
124075396ef9SLisandro Dalcin     ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr);
124145fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
124245fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
124345fc7adcSBarry Smith   }
124432a4b47aSMatthew Knepley #endif
124545fc7adcSBarry Smith 
12467adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1247c1f60f51SBarry Smith   /*
1248dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1249c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1250c1f60f51SBarry Smith    */
125131615d04SLois Curfman McInnes   if (flg) {
1252272ac6f2SLois Curfman McInnes     Mat  J;
1253b5d62d44SBarry Smith     KSP ksp;
125494b7f48cSBarry Smith     PC   pc;
125575396ef9SLisandro Dalcin     /* create and set matrix-free operator */
1256fef1beadSBarry Smith     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
12573ec795f1SBarry Smith     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
125875396ef9SLisandro Dalcin     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
12593ec795f1SBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
12603a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1261f3ef73ceSBarry Smith     /* force no preconditioner */
126294b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1263b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1264a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1265a9815358SBarry Smith     if (!flg) {
126675396ef9SLisandro Dalcin       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
1267f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1268272ac6f2SLois Curfman McInnes     }
1269a9815358SBarry Smith   }
1270f3ef73ceSBarry Smith 
127185385478SLisandro Dalcin   if (!snes->vec_func && !snes->vec_rhs) {
12721096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
12731096aae1SMatthew Knepley   }
127485385478SLisandro Dalcin   if (!snes->ops->computefunction && !snes->vec_rhs) {
12751096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
12761096aae1SMatthew Knepley   }
127775396ef9SLisandro Dalcin   if (!snes->jacobian) {
127875396ef9SLisandro Dalcin     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
127975396ef9SLisandro Dalcin   }
1280a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
128129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1282a8c6a408SBarry Smith   }
1283a703fe33SLois Curfman McInnes 
1284410397dcSLisandro Dalcin   if (snes->ops->setup) {
1285410397dcSLisandro Dalcin     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1286410397dcSLisandro Dalcin   }
12877aaed0d8SBarry Smith   snes->setupcalled = PETSC_TRUE;
12883a40ed3dSBarry Smith   PetscFunctionReturn(0);
12899b94acceSBarry Smith }
12909b94acceSBarry Smith 
12914a2ae208SSatish Balay #undef __FUNCT__
12924a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
129352baeb72SSatish Balay /*@
12949b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
12959b94acceSBarry Smith    with SNESCreate().
12969b94acceSBarry Smith 
1297c7afd0dbSLois Curfman McInnes    Collective on SNES
1298c7afd0dbSLois Curfman McInnes 
12999b94acceSBarry Smith    Input Parameter:
13009b94acceSBarry Smith .  snes - the SNES context
13019b94acceSBarry Smith 
130236851e7fSLois Curfman McInnes    Level: beginner
130336851e7fSLois Curfman McInnes 
13049b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
13059b94acceSBarry Smith 
130663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
13079b94acceSBarry Smith @*/
130863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
13099b94acceSBarry Smith {
13106849ba73SBarry Smith   PetscErrorCode ierr;
13113a40ed3dSBarry Smith 
13123a40ed3dSBarry Smith   PetscFunctionBegin;
13134482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13147adad957SLisandro Dalcin   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1315d4bb536fSBarry Smith 
1316be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
13170f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1318e7788613SBarry Smith   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
131985385478SLisandro Dalcin 
132085385478SLisandro Dalcin   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
132185385478SLisandro Dalcin   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
132285385478SLisandro Dalcin   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
13233a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
13243a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
132594b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
132685385478SLisandro Dalcin   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1327522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1328a6570f20SBarry Smith   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
13297f7931b9SBarry Smith   if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);}
1330a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
13313a40ed3dSBarry Smith   PetscFunctionReturn(0);
13329b94acceSBarry Smith }
13339b94acceSBarry Smith 
13349b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
13359b94acceSBarry Smith 
13364a2ae208SSatish Balay #undef __FUNCT__
1337*a8054027SBarry Smith #define __FUNCT__ "SNESSetLagPreconditioner"
1338*a8054027SBarry Smith /*@
1339*a8054027SBarry Smith    SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve.
1340*a8054027SBarry Smith 
1341*a8054027SBarry Smith    Collective on SNES
1342*a8054027SBarry Smith 
1343*a8054027SBarry Smith    Input Parameters:
1344*a8054027SBarry Smith +  snes - the SNES context
1345*a8054027SBarry 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
1346*a8054027SBarry Smith          the Jacobian is built etc.
1347*a8054027SBarry Smith 
1348*a8054027SBarry Smith    Options Database Keys:
1349*a8054027SBarry Smith .    -snes_lag_preconditioner <lag>
1350*a8054027SBarry Smith 
1351*a8054027SBarry Smith    Notes:
1352*a8054027SBarry Smith    The default is 1
1353*a8054027SBarry Smith    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1354*a8054027SBarry Smith    If  -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use
1355*a8054027SBarry Smith 
1356*a8054027SBarry Smith    Level: intermediate
1357*a8054027SBarry Smith 
1358*a8054027SBarry Smith .keywords: SNES, nonlinear, set, convergence, tolerances
1359*a8054027SBarry Smith 
1360*a8054027SBarry Smith .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner()
1361*a8054027SBarry Smith 
1362*a8054027SBarry Smith @*/
1363*a8054027SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag)
1364*a8054027SBarry Smith {
1365*a8054027SBarry Smith   PetscFunctionBegin;
1366*a8054027SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1367*a8054027SBarry Smith   if (lag < -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -1, 1 or greater");
1368*a8054027SBarry Smith   if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0");
1369*a8054027SBarry Smith   snes->lagpreconditioner = lag;
1370*a8054027SBarry Smith   PetscFunctionReturn(0);
1371*a8054027SBarry Smith }
1372*a8054027SBarry Smith 
1373*a8054027SBarry Smith #undef __FUNCT__
1374*a8054027SBarry Smith #define __FUNCT__ "SNESGetLagPreconditioner"
1375*a8054027SBarry Smith /*@
1376*a8054027SBarry Smith    SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt
1377*a8054027SBarry Smith 
1378*a8054027SBarry Smith    Collective on SNES
1379*a8054027SBarry Smith 
1380*a8054027SBarry Smith    Input Parameter:
1381*a8054027SBarry Smith .  snes - the SNES context
1382*a8054027SBarry Smith 
1383*a8054027SBarry Smith    Output Parameter:
1384*a8054027SBarry 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
1385*a8054027SBarry Smith          the Jacobian is built etc.
1386*a8054027SBarry Smith 
1387*a8054027SBarry Smith    Options Database Keys:
1388*a8054027SBarry Smith .    -snes_lag_preconditioner <lag>
1389*a8054027SBarry Smith 
1390*a8054027SBarry Smith    Notes:
1391*a8054027SBarry Smith    The default is 1
1392*a8054027SBarry Smith    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1
1393*a8054027SBarry Smith 
1394*a8054027SBarry Smith    Level: intermediate
1395*a8054027SBarry Smith 
1396*a8054027SBarry Smith .keywords: SNES, nonlinear, set, convergence, tolerances
1397*a8054027SBarry Smith 
1398*a8054027SBarry Smith .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner()
1399*a8054027SBarry Smith 
1400*a8054027SBarry Smith @*/
1401*a8054027SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag)
1402*a8054027SBarry Smith {
1403*a8054027SBarry Smith   PetscFunctionBegin;
1404*a8054027SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1405*a8054027SBarry Smith   *lag = snes->lagpreconditioner;
1406*a8054027SBarry Smith   PetscFunctionReturn(0);
1407*a8054027SBarry Smith }
1408*a8054027SBarry Smith 
1409*a8054027SBarry Smith #undef __FUNCT__
14104a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
14119b94acceSBarry Smith /*@
1412d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
14139b94acceSBarry Smith 
1414c7afd0dbSLois Curfman McInnes    Collective on SNES
1415c7afd0dbSLois Curfman McInnes 
14169b94acceSBarry Smith    Input Parameters:
1417c7afd0dbSLois Curfman McInnes +  snes - the SNES context
141870441072SBarry Smith .  abstol - absolute convergence tolerance
141933174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
142033174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
142133174efeSLois Curfman McInnes            of the change in the solution between steps
142233174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1423c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1424fee21e36SBarry Smith 
142533174efeSLois Curfman McInnes    Options Database Keys:
142670441072SBarry Smith +    -snes_atol <abstol> - Sets abstol
1427c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1428c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1429c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1430c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
14319b94acceSBarry Smith 
1432d7a720efSLois Curfman McInnes    Notes:
14339b94acceSBarry Smith    The default maximum number of iterations is 50.
14349b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
14359b94acceSBarry Smith 
143636851e7fSLois Curfman McInnes    Level: intermediate
143736851e7fSLois Curfman McInnes 
143833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
14399b94acceSBarry Smith 
14402492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance()
14419b94acceSBarry Smith @*/
144263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
14439b94acceSBarry Smith {
14443a40ed3dSBarry Smith   PetscFunctionBegin;
14454482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
144670441072SBarry Smith   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1447d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1448d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1449d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1450d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
14513a40ed3dSBarry Smith   PetscFunctionReturn(0);
14529b94acceSBarry Smith }
14539b94acceSBarry Smith 
14544a2ae208SSatish Balay #undef __FUNCT__
14554a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
14569b94acceSBarry Smith /*@
145733174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
145833174efeSLois Curfman McInnes 
1459c7afd0dbSLois Curfman McInnes    Not Collective
1460c7afd0dbSLois Curfman McInnes 
146133174efeSLois Curfman McInnes    Input Parameters:
1462c7afd0dbSLois Curfman McInnes +  snes - the SNES context
146385385478SLisandro Dalcin .  atol - absolute convergence tolerance
146433174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
146533174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
146633174efeSLois Curfman McInnes            of the change in the solution between steps
146733174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1468c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1469fee21e36SBarry Smith 
147033174efeSLois Curfman McInnes    Notes:
147133174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
147233174efeSLois Curfman McInnes 
147336851e7fSLois Curfman McInnes    Level: intermediate
147436851e7fSLois Curfman McInnes 
147533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
147633174efeSLois Curfman McInnes 
147733174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
147833174efeSLois Curfman McInnes @*/
147985385478SLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
148033174efeSLois Curfman McInnes {
14813a40ed3dSBarry Smith   PetscFunctionBegin;
14824482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
148385385478SLisandro Dalcin   if (atol)  *atol  = snes->abstol;
148433174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
148533174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
148633174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
148733174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
14883a40ed3dSBarry Smith   PetscFunctionReturn(0);
148933174efeSLois Curfman McInnes }
149033174efeSLois Curfman McInnes 
14914a2ae208SSatish Balay #undef __FUNCT__
14924a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
149333174efeSLois Curfman McInnes /*@
14949b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
14959b94acceSBarry Smith 
1496fee21e36SBarry Smith    Collective on SNES
1497fee21e36SBarry Smith 
1498c7afd0dbSLois Curfman McInnes    Input Parameters:
1499c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1500c7afd0dbSLois Curfman McInnes -  tol - tolerance
1501c7afd0dbSLois Curfman McInnes 
15029b94acceSBarry Smith    Options Database Key:
1503c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
15049b94acceSBarry Smith 
150536851e7fSLois Curfman McInnes    Level: intermediate
150636851e7fSLois Curfman McInnes 
15079b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
15089b94acceSBarry Smith 
15092492ecdbSBarry Smith .seealso: SNESSetTolerances()
15109b94acceSBarry Smith @*/
151163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
15129b94acceSBarry Smith {
15133a40ed3dSBarry Smith   PetscFunctionBegin;
15144482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
15159b94acceSBarry Smith   snes->deltatol = tol;
15163a40ed3dSBarry Smith   PetscFunctionReturn(0);
15179b94acceSBarry Smith }
15189b94acceSBarry Smith 
1519df9fa365SBarry Smith /*
1520df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1521df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1522df9fa365SBarry Smith    macros instead of functions
1523df9fa365SBarry Smith */
15244a2ae208SSatish Balay #undef __FUNCT__
1525a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLG"
1526a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1527ce1608b8SBarry Smith {
1528dfbe8321SBarry Smith   PetscErrorCode ierr;
1529ce1608b8SBarry Smith 
1530ce1608b8SBarry Smith   PetscFunctionBegin;
15314482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1532a6570f20SBarry Smith   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1533ce1608b8SBarry Smith   PetscFunctionReturn(0);
1534ce1608b8SBarry Smith }
1535ce1608b8SBarry Smith 
15364a2ae208SSatish Balay #undef __FUNCT__
1537a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLGCreate"
1538a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1539df9fa365SBarry Smith {
1540dfbe8321SBarry Smith   PetscErrorCode ierr;
1541df9fa365SBarry Smith 
1542df9fa365SBarry Smith   PetscFunctionBegin;
1543a6570f20SBarry Smith   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1544df9fa365SBarry Smith   PetscFunctionReturn(0);
1545df9fa365SBarry Smith }
1546df9fa365SBarry Smith 
15474a2ae208SSatish Balay #undef __FUNCT__
1548a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLGDestroy"
1549a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1550df9fa365SBarry Smith {
1551dfbe8321SBarry Smith   PetscErrorCode ierr;
1552df9fa365SBarry Smith 
1553df9fa365SBarry Smith   PetscFunctionBegin;
1554a6570f20SBarry Smith   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1555df9fa365SBarry Smith   PetscFunctionReturn(0);
1556df9fa365SBarry Smith }
1557df9fa365SBarry Smith 
15589b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
15599b94acceSBarry Smith 
15604a2ae208SSatish Balay #undef __FUNCT__
1561a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorSet"
15629b94acceSBarry Smith /*@C
1563a6570f20SBarry Smith    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
15649b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
15659b94acceSBarry Smith    progress.
15669b94acceSBarry Smith 
1567fee21e36SBarry Smith    Collective on SNES
1568fee21e36SBarry Smith 
1569c7afd0dbSLois Curfman McInnes    Input Parameters:
1570c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1571c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1572b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1573e8105e01SRichard Katz           monitor routine (use PETSC_NULL if no context is desired)
1574b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1575b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
15769b94acceSBarry Smith 
1577c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1578a7cc72afSBarry Smith $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1579c7afd0dbSLois Curfman McInnes 
1580c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1581c7afd0dbSLois Curfman McInnes .    its - iteration number
1582c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
158340a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
15849b94acceSBarry Smith 
15859665c990SLois Curfman McInnes    Options Database Keys:
1586a6570f20SBarry Smith +    -snes_monitor        - sets SNESMonitorDefault()
1587a6570f20SBarry Smith .    -snes_monitor_draw    - sets line graph monitor,
1588a6570f20SBarry Smith                             uses SNESMonitorLGCreate()
1589a6570f20SBarry Smith _    -snes_monitor_cancel - cancels all monitors that have
1590c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1591a6570f20SBarry Smith                             calls to SNESMonitorSet(), but
1592c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1593c7afd0dbSLois Curfman McInnes                             the options database.
15949665c990SLois Curfman McInnes 
1595639f9d9dSBarry Smith    Notes:
15966bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
1597a6570f20SBarry Smith    SNESMonitorSet() multiple times; all will be called in the
15986bc08f3fSLois Curfman McInnes    order in which they were set.
1599639f9d9dSBarry Smith 
1600025f1a04SBarry Smith    Fortran notes: Only a single monitor function can be set for each SNES object
1601025f1a04SBarry Smith 
160236851e7fSLois Curfman McInnes    Level: intermediate
160336851e7fSLois Curfman McInnes 
16049b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
16059b94acceSBarry Smith 
1606a6570f20SBarry Smith .seealso: SNESMonitorDefault(), SNESMonitorCancel()
16079b94acceSBarry Smith @*/
1608b90d0a6eSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
16099b94acceSBarry Smith {
1610b90d0a6eSBarry Smith   PetscInt i;
1611b90d0a6eSBarry Smith 
16123a40ed3dSBarry Smith   PetscFunctionBegin;
16134482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1614639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
161529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1616639f9d9dSBarry Smith   }
1617b90d0a6eSBarry Smith   for (i=0; i<snes->numbermonitors;i++) {
1618b90d0a6eSBarry Smith     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1619b90d0a6eSBarry Smith 
1620b90d0a6eSBarry Smith     /* check if both default monitors that share common ASCII viewer */
1621b90d0a6eSBarry Smith     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1622b90d0a6eSBarry Smith       if (mctx && snes->monitorcontext[i]) {
1623b90d0a6eSBarry Smith         PetscErrorCode          ierr;
1624b90d0a6eSBarry Smith         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1625b90d0a6eSBarry Smith         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1626b90d0a6eSBarry Smith         if (viewer1->viewer == viewer2->viewer) {
1627b90d0a6eSBarry Smith           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1628b90d0a6eSBarry Smith           PetscFunctionReturn(0);
1629b90d0a6eSBarry Smith         }
1630b90d0a6eSBarry Smith       }
1631b90d0a6eSBarry Smith     }
1632b90d0a6eSBarry Smith   }
1633b90d0a6eSBarry Smith   snes->monitor[snes->numbermonitors]           = monitor;
1634b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1635639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
16363a40ed3dSBarry Smith   PetscFunctionReturn(0);
16379b94acceSBarry Smith }
16389b94acceSBarry Smith 
16394a2ae208SSatish Balay #undef __FUNCT__
1640a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorCancel"
16415cd90555SBarry Smith /*@C
1642a6570f20SBarry Smith    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
16435cd90555SBarry Smith 
1644c7afd0dbSLois Curfman McInnes    Collective on SNES
1645c7afd0dbSLois Curfman McInnes 
16465cd90555SBarry Smith    Input Parameters:
16475cd90555SBarry Smith .  snes - the SNES context
16485cd90555SBarry Smith 
16491a480d89SAdministrator    Options Database Key:
1650a6570f20SBarry Smith .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1651a6570f20SBarry Smith     into a code by calls to SNESMonitorSet(), but does not cancel those
1652c7afd0dbSLois Curfman McInnes     set via the options database
16535cd90555SBarry Smith 
16545cd90555SBarry Smith    Notes:
16555cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
16565cd90555SBarry Smith 
165736851e7fSLois Curfman McInnes    Level: intermediate
165836851e7fSLois Curfman McInnes 
16595cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
16605cd90555SBarry Smith 
1661a6570f20SBarry Smith .seealso: SNESMonitorDefault(), SNESMonitorSet()
16625cd90555SBarry Smith @*/
1663a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
16645cd90555SBarry Smith {
1665d952e501SBarry Smith   PetscErrorCode ierr;
1666d952e501SBarry Smith   PetscInt       i;
1667d952e501SBarry Smith 
16685cd90555SBarry Smith   PetscFunctionBegin;
16694482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1670d952e501SBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1671d952e501SBarry Smith     if (snes->monitordestroy[i]) {
1672d952e501SBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1673d952e501SBarry Smith     }
1674d952e501SBarry Smith   }
16755cd90555SBarry Smith   snes->numbermonitors = 0;
16765cd90555SBarry Smith   PetscFunctionReturn(0);
16775cd90555SBarry Smith }
16785cd90555SBarry Smith 
16794a2ae208SSatish Balay #undef __FUNCT__
16804a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
16819b94acceSBarry Smith /*@C
16829b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
16839b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
16849b94acceSBarry Smith 
1685fee21e36SBarry Smith    Collective on SNES
1686fee21e36SBarry Smith 
1687c7afd0dbSLois Curfman McInnes    Input Parameters:
1688c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1689c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
16907f7931b9SBarry Smith .  cctx - [optional] context for private data for the convergence routine  (may be PETSC_NULL)
16917f7931b9SBarry Smith -  destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran)
16929b94acceSBarry Smith 
1693c7afd0dbSLois Curfman McInnes    Calling sequence of func:
169406ee9f85SBarry Smith $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1695c7afd0dbSLois Curfman McInnes 
1696c7afd0dbSLois Curfman McInnes +    snes - the SNES context
169706ee9f85SBarry Smith .    it - current iteration (0 is the first and is before any Newton step)
1698c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1699184914b5SBarry Smith .    reason - reason for convergence/divergence
1700c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
17014b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
17024b27c08aSLois Curfman McInnes -    f - 2-norm of function
17039b94acceSBarry Smith 
170436851e7fSLois Curfman McInnes    Level: advanced
170536851e7fSLois Curfman McInnes 
17069b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
17079b94acceSBarry Smith 
170885385478SLisandro Dalcin .seealso: SNESDefaultConverged(), SNESSkipConverged()
17099b94acceSBarry Smith @*/
17107f7931b9SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
17119b94acceSBarry Smith {
17127f7931b9SBarry Smith   PetscErrorCode ierr;
17137f7931b9SBarry Smith 
17143a40ed3dSBarry Smith   PetscFunctionBegin;
17154482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
171685385478SLisandro Dalcin   if (!func) func = SNESSkipConverged;
17177f7931b9SBarry Smith   if (snes->ops->convergeddestroy) {
17187f7931b9SBarry Smith     ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);
17197f7931b9SBarry Smith   }
172085385478SLisandro Dalcin   snes->ops->converged        = func;
17217f7931b9SBarry Smith   snes->ops->convergeddestroy = destroy;
172285385478SLisandro Dalcin   snes->cnvP                  = cctx;
17233a40ed3dSBarry Smith   PetscFunctionReturn(0);
17249b94acceSBarry Smith }
17259b94acceSBarry Smith 
17264a2ae208SSatish Balay #undef __FUNCT__
17274a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
172852baeb72SSatish Balay /*@
1729184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1730184914b5SBarry Smith 
1731184914b5SBarry Smith    Not Collective
1732184914b5SBarry Smith 
1733184914b5SBarry Smith    Input Parameter:
1734184914b5SBarry Smith .  snes - the SNES context
1735184914b5SBarry Smith 
1736184914b5SBarry Smith    Output Parameter:
1737e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1738184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1739184914b5SBarry Smith 
1740184914b5SBarry Smith    Level: intermediate
1741184914b5SBarry Smith 
1742184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1743184914b5SBarry Smith 
1744184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1745184914b5SBarry Smith 
174685385478SLisandro Dalcin .seealso: SNESSetConvergenceTest(), SNESConvergedReason
1747184914b5SBarry Smith @*/
174863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1749184914b5SBarry Smith {
1750184914b5SBarry Smith   PetscFunctionBegin;
17514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17524482741eSBarry Smith   PetscValidPointer(reason,2);
1753184914b5SBarry Smith   *reason = snes->reason;
1754184914b5SBarry Smith   PetscFunctionReturn(0);
1755184914b5SBarry Smith }
1756184914b5SBarry Smith 
17574a2ae208SSatish Balay #undef __FUNCT__
17584a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1759c9005455SLois Curfman McInnes /*@
1760c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1761c9005455SLois Curfman McInnes 
1762fee21e36SBarry Smith    Collective on SNES
1763fee21e36SBarry Smith 
1764c7afd0dbSLois Curfman McInnes    Input Parameters:
1765c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1766c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1767cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1768758f92a0SBarry Smith .  na  - size of a and its
176964731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1770758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1771c7afd0dbSLois Curfman McInnes 
1772c9005455SLois Curfman McInnes    Notes:
17734b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1774c9005455SLois Curfman McInnes    at each step.
1775c9005455SLois Curfman McInnes 
1776c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1777c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1778c9005455SLois Curfman McInnes    during the section of code that is being timed.
1779c9005455SLois Curfman McInnes 
178036851e7fSLois Curfman McInnes    Level: intermediate
178136851e7fSLois Curfman McInnes 
1782c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1783758f92a0SBarry Smith 
178408405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1785758f92a0SBarry Smith 
1786c9005455SLois Curfman McInnes @*/
1787a562a398SLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset)
1788c9005455SLois Curfman McInnes {
17893a40ed3dSBarry Smith   PetscFunctionBegin;
17904482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17914482741eSBarry Smith   if (na)  PetscValidScalarPointer(a,2);
1792a562a398SLisandro Dalcin   if (its) PetscValidIntPointer(its,3);
1793c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1794758f92a0SBarry Smith   snes->conv_hist_its   = its;
1795758f92a0SBarry Smith   snes->conv_hist_max   = na;
1796a12bc48fSLisandro Dalcin   snes->conv_hist_len   = 0;
1797758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1798758f92a0SBarry Smith   PetscFunctionReturn(0);
1799758f92a0SBarry Smith }
1800758f92a0SBarry Smith 
18014a2ae208SSatish Balay #undef __FUNCT__
18024a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
18030c4c9dddSBarry Smith /*@C
1804758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1805758f92a0SBarry Smith 
1806758f92a0SBarry Smith    Collective on SNES
1807758f92a0SBarry Smith 
1808758f92a0SBarry Smith    Input Parameter:
1809758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1810758f92a0SBarry Smith 
1811758f92a0SBarry Smith    Output Parameters:
1812758f92a0SBarry Smith .  a   - array to hold history
1813758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1814758f92a0SBarry Smith          negative if not converged) for each solve.
1815758f92a0SBarry Smith -  na  - size of a and its
1816758f92a0SBarry Smith 
1817758f92a0SBarry Smith    Notes:
1818758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1819758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1820758f92a0SBarry Smith 
1821758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1822758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1823758f92a0SBarry Smith    during the section of code that is being timed.
1824758f92a0SBarry Smith 
1825758f92a0SBarry Smith    Level: intermediate
1826758f92a0SBarry Smith 
1827758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1828758f92a0SBarry Smith 
1829758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1830758f92a0SBarry Smith 
1831758f92a0SBarry Smith @*/
183263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1833758f92a0SBarry Smith {
1834758f92a0SBarry Smith   PetscFunctionBegin;
18354482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1836758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1837758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1838758f92a0SBarry Smith   if (na)  *na  = snes->conv_hist_len;
18393a40ed3dSBarry Smith   PetscFunctionReturn(0);
1840c9005455SLois Curfman McInnes }
1841c9005455SLois Curfman McInnes 
1842e74ef692SMatthew Knepley #undef __FUNCT__
1843e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1844ac226902SBarry Smith /*@C
184576b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
18467e4bb74cSBarry Smith   at the beginning o every iteration of the nonlinear solve. Specifically
18477e4bb74cSBarry Smith   it is called just before the Jacobian is "evaluated".
184876b2cf59SMatthew Knepley 
184976b2cf59SMatthew Knepley   Collective on SNES
185076b2cf59SMatthew Knepley 
185176b2cf59SMatthew Knepley   Input Parameters:
185276b2cf59SMatthew Knepley . snes - The nonlinear solver context
185376b2cf59SMatthew Knepley . func - The function
185476b2cf59SMatthew Knepley 
185576b2cf59SMatthew Knepley   Calling sequence of func:
1856b5d30489SBarry Smith . func (SNES snes, PetscInt step);
185776b2cf59SMatthew Knepley 
185876b2cf59SMatthew Knepley . step - The current step of the iteration
185976b2cf59SMatthew Knepley 
186076b2cf59SMatthew Knepley   Level: intermediate
186176b2cf59SMatthew Knepley 
186276b2cf59SMatthew Knepley .keywords: SNES, update
1863b5d30489SBarry Smith 
186485385478SLisandro Dalcin .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
186576b2cf59SMatthew Knepley @*/
186663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
186776b2cf59SMatthew Knepley {
186876b2cf59SMatthew Knepley   PetscFunctionBegin;
18694482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
1870e7788613SBarry Smith   snes->ops->update = func;
187176b2cf59SMatthew Knepley   PetscFunctionReturn(0);
187276b2cf59SMatthew Knepley }
187376b2cf59SMatthew Knepley 
1874e74ef692SMatthew Knepley #undef __FUNCT__
1875e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
187676b2cf59SMatthew Knepley /*@
187776b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
187876b2cf59SMatthew Knepley 
187976b2cf59SMatthew Knepley   Not collective
188076b2cf59SMatthew Knepley 
188176b2cf59SMatthew Knepley   Input Parameters:
188276b2cf59SMatthew Knepley . snes - The nonlinear solver context
188376b2cf59SMatthew Knepley . step - The current step of the iteration
188476b2cf59SMatthew Knepley 
1885205452f4SMatthew Knepley   Level: intermediate
1886205452f4SMatthew Knepley 
188776b2cf59SMatthew Knepley .keywords: SNES, update
1888a6570f20SBarry Smith .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
188976b2cf59SMatthew Knepley @*/
189063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
189176b2cf59SMatthew Knepley {
189276b2cf59SMatthew Knepley   PetscFunctionBegin;
189376b2cf59SMatthew Knepley   PetscFunctionReturn(0);
189476b2cf59SMatthew Knepley }
189576b2cf59SMatthew Knepley 
18964a2ae208SSatish Balay #undef __FUNCT__
18974a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
18989b94acceSBarry Smith /*
18999b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
19009b94acceSBarry Smith    positive parameter delta.
19019b94acceSBarry Smith 
19029b94acceSBarry Smith     Input Parameters:
1903c7afd0dbSLois Curfman McInnes +   snes - the SNES context
19049b94acceSBarry Smith .   y - approximate solution of linear system
19059b94acceSBarry Smith .   fnorm - 2-norm of current function
1906c7afd0dbSLois Curfman McInnes -   delta - trust region size
19079b94acceSBarry Smith 
19089b94acceSBarry Smith     Output Parameters:
1909c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
19109b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
19119b94acceSBarry Smith     region, and exceeds zero otherwise.
1912c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
19139b94acceSBarry Smith 
19149b94acceSBarry Smith     Note:
19154b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
19169b94acceSBarry Smith     is set to be the maximum allowable step size.
19179b94acceSBarry Smith 
19189b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
19199b94acceSBarry Smith */
1920dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
19219b94acceSBarry Smith {
1922064f8208SBarry Smith   PetscReal      nrm;
1923ea709b57SSatish Balay   PetscScalar    cnorm;
1924dfbe8321SBarry Smith   PetscErrorCode ierr;
19253a40ed3dSBarry Smith 
19263a40ed3dSBarry Smith   PetscFunctionBegin;
19274482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19284482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1929c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1930184914b5SBarry Smith 
1931064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1932064f8208SBarry Smith   if (nrm > *delta) {
1933064f8208SBarry Smith      nrm = *delta/nrm;
1934064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1935064f8208SBarry Smith      cnorm = nrm;
19362dcb1b2aSMatthew Knepley      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
19379b94acceSBarry Smith      *ynorm = *delta;
19389b94acceSBarry Smith   } else {
19399b94acceSBarry Smith      *gpnorm = 0.0;
1940064f8208SBarry Smith      *ynorm = nrm;
19419b94acceSBarry Smith   }
19423a40ed3dSBarry Smith   PetscFunctionReturn(0);
19439b94acceSBarry Smith }
19449b94acceSBarry Smith 
19454a2ae208SSatish Balay #undef __FUNCT__
19464a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
19476ce558aeSBarry Smith /*@C
1948f69a0ea3SMatthew Knepley    SNESSolve - Solves a nonlinear system F(x) = b.
1949f69a0ea3SMatthew Knepley    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
19509b94acceSBarry Smith 
1951c7afd0dbSLois Curfman McInnes    Collective on SNES
1952c7afd0dbSLois Curfman McInnes 
1953b2002411SLois Curfman McInnes    Input Parameters:
1954c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1955f69a0ea3SMatthew Knepley .  b - the constant part of the equation, or PETSC_NULL to use zero.
195685385478SLisandro Dalcin -  x - the solution vector.
19579b94acceSBarry Smith 
1958b2002411SLois Curfman McInnes    Notes:
19598ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
19608ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
19618ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
19628ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
19638ddd3da0SLois Curfman McInnes 
196436851e7fSLois Curfman McInnes    Level: beginner
196536851e7fSLois Curfman McInnes 
19669b94acceSBarry Smith .keywords: SNES, nonlinear, solve
19679b94acceSBarry Smith 
196885385478SLisandro Dalcin .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
19699b94acceSBarry Smith @*/
1970f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
19719b94acceSBarry Smith {
1972dfbe8321SBarry Smith   PetscErrorCode ierr;
1973f1af5d2fSBarry Smith   PetscTruth     flg;
1974eabae89aSBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
1975eabae89aSBarry Smith   PetscViewer    viewer;
1976052efed2SBarry Smith 
19773a40ed3dSBarry Smith   PetscFunctionBegin;
19784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1979f69a0ea3SMatthew Knepley   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
1980f69a0ea3SMatthew Knepley   PetscCheckSameComm(snes,1,x,3);
198185385478SLisandro Dalcin   if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2);
198285385478SLisandro Dalcin   if (b) PetscCheckSameComm(snes,1,b,2);
198385385478SLisandro Dalcin 
198485385478SLisandro Dalcin   /* set solution vector */
198585385478SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
198685385478SLisandro Dalcin   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
198785385478SLisandro Dalcin   snes->vec_sol = x;
198885385478SLisandro Dalcin   /* set afine vector if provided */
198985385478SLisandro Dalcin   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
199085385478SLisandro Dalcin   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
199185385478SLisandro Dalcin   snes->vec_rhs = b;
199285385478SLisandro Dalcin 
199385385478SLisandro Dalcin   if (!snes->vec_func && snes->vec_rhs) {
199485385478SLisandro Dalcin     ierr = VecDuplicate(b, &snes->vec_func); CHKERRQ(ierr);
199570e92668SMatthew Knepley   }
19963f149594SLisandro Dalcin 
199770e92668SMatthew Knepley   ierr = SNESSetUp(snes);CHKERRQ(ierr);
19983f149594SLisandro Dalcin 
1999abc0a331SBarry Smith   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
200050ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
2001d5e45103SBarry Smith 
20023f149594SLisandro Dalcin   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
20034936397dSBarry Smith   ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr);
200485385478SLisandro Dalcin   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
20054936397dSBarry Smith   if (snes->domainerror){
20064936397dSBarry Smith     snes->reason      = SNES_DIVERGED_FUNCTION_DOMAIN;
20074936397dSBarry Smith     snes->domainerror = PETSC_FALSE;
20084936397dSBarry Smith   }
200985385478SLisandro Dalcin 
20103f149594SLisandro Dalcin   if (!snes->reason) {
20113f149594SLisandro Dalcin     SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
20123f149594SLisandro Dalcin   }
20133f149594SLisandro Dalcin 
20147adad957SLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
2015eabae89aSBarry Smith   if (flg && !PetscPreLoadingOn) {
20167adad957SLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
2017eabae89aSBarry Smith     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
2018eabae89aSBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
2019eabae89aSBarry Smith   }
2020eabae89aSBarry Smith 
20217adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
2022da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
20235968eb51SBarry Smith   if (snes->printreason) {
20245968eb51SBarry Smith     if (snes->reason > 0) {
20257adad957SLisandro Dalcin       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
20265968eb51SBarry Smith     } else {
20277adad957SLisandro Dalcin       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
20285968eb51SBarry Smith     }
20295968eb51SBarry Smith   }
20305968eb51SBarry Smith 
20313a40ed3dSBarry Smith   PetscFunctionReturn(0);
20329b94acceSBarry Smith }
20339b94acceSBarry Smith 
20349b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
20359b94acceSBarry Smith 
20364a2ae208SSatish Balay #undef __FUNCT__
20374a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
203882bf6240SBarry Smith /*@C
20394b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
20409b94acceSBarry Smith 
2041fee21e36SBarry Smith    Collective on SNES
2042fee21e36SBarry Smith 
2043c7afd0dbSLois Curfman McInnes    Input Parameters:
2044c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2045454a90a3SBarry Smith -  type - a known method
2046c7afd0dbSLois Curfman McInnes 
2047c7afd0dbSLois Curfman McInnes    Options Database Key:
2048454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
2049c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
2050ae12b187SLois Curfman McInnes 
20519b94acceSBarry Smith    Notes:
2052e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
20534b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
2054c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
20554b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
2056c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
20579b94acceSBarry Smith 
2058ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
2059ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
2060ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
2061ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
2062ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
2063ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
2064ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
2065ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
2066ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
2067b0a32e0cSBarry Smith   appropriate method.
206836851e7fSLois Curfman McInnes 
206936851e7fSLois Curfman McInnes   Level: intermediate
2070a703fe33SLois Curfman McInnes 
2071454a90a3SBarry Smith .keywords: SNES, set, type
2072435da068SBarry Smith 
2073435da068SBarry Smith .seealso: SNESType, SNESCreate()
2074435da068SBarry Smith 
20759b94acceSBarry Smith @*/
2076a313700dSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type)
20779b94acceSBarry Smith {
2078dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
20796831982aSBarry Smith   PetscTruth     match;
20803a40ed3dSBarry Smith 
20813a40ed3dSBarry Smith   PetscFunctionBegin;
20824482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
20834482741eSBarry Smith   PetscValidCharPointer(type,2);
208482bf6240SBarry Smith 
20856831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
20860f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
208792ff6ae8SBarry Smith 
20887adad957SLisandro Dalcin   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
2089958c9bccSBarry Smith   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
209075396ef9SLisandro Dalcin   /* Destroy the previous private SNES context */
209175396ef9SLisandro Dalcin   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
209275396ef9SLisandro Dalcin   /* Reinitialize function pointers in SNESOps structure */
209375396ef9SLisandro Dalcin   snes->ops->setup          = 0;
209475396ef9SLisandro Dalcin   snes->ops->solve          = 0;
209575396ef9SLisandro Dalcin   snes->ops->view           = 0;
209675396ef9SLisandro Dalcin   snes->ops->setfromoptions = 0;
209775396ef9SLisandro Dalcin   snes->ops->destroy        = 0;
209875396ef9SLisandro Dalcin   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
209975396ef9SLisandro Dalcin   snes->setupcalled = PETSC_FALSE;
21003a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
2101454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
21023a40ed3dSBarry Smith   PetscFunctionReturn(0);
21039b94acceSBarry Smith }
21049b94acceSBarry Smith 
2105a847f771SSatish Balay 
21069b94acceSBarry Smith /* --------------------------------------------------------------------- */
21074a2ae208SSatish Balay #undef __FUNCT__
21084a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
210952baeb72SSatish Balay /*@
21109b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2111f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
21129b94acceSBarry Smith 
2113fee21e36SBarry Smith    Not Collective
2114fee21e36SBarry Smith 
211536851e7fSLois Curfman McInnes    Level: advanced
211636851e7fSLois Curfman McInnes 
21179b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
21189b94acceSBarry Smith 
21199b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
21209b94acceSBarry Smith @*/
212163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
21229b94acceSBarry Smith {
2123dfbe8321SBarry Smith   PetscErrorCode ierr;
212482bf6240SBarry Smith 
21253a40ed3dSBarry Smith   PetscFunctionBegin;
21261441b1d3SBarry Smith   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
21274c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
21283a40ed3dSBarry Smith   PetscFunctionReturn(0);
21299b94acceSBarry Smith }
21309b94acceSBarry Smith 
21314a2ae208SSatish Balay #undef __FUNCT__
21324a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
21339b94acceSBarry Smith /*@C
21349a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
21359b94acceSBarry Smith 
2136c7afd0dbSLois Curfman McInnes    Not Collective
2137c7afd0dbSLois Curfman McInnes 
21389b94acceSBarry Smith    Input Parameter:
21394b0e389bSBarry Smith .  snes - nonlinear solver context
21409b94acceSBarry Smith 
21419b94acceSBarry Smith    Output Parameter:
21423a7fca6bSBarry Smith .  type - SNES method (a character string)
21439b94acceSBarry Smith 
214436851e7fSLois Curfman McInnes    Level: intermediate
214536851e7fSLois Curfman McInnes 
2146454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
21479b94acceSBarry Smith @*/
2148a313700dSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type)
21499b94acceSBarry Smith {
21503a40ed3dSBarry Smith   PetscFunctionBegin;
21514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
21524482741eSBarry Smith   PetscValidPointer(type,2);
21537adad957SLisandro Dalcin   *type = ((PetscObject)snes)->type_name;
21543a40ed3dSBarry Smith   PetscFunctionReturn(0);
21559b94acceSBarry Smith }
21569b94acceSBarry Smith 
21574a2ae208SSatish Balay #undef __FUNCT__
21584a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
215952baeb72SSatish Balay /*@
21609b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
21619b94acceSBarry Smith    stored.
21629b94acceSBarry Smith 
2163c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2164c7afd0dbSLois Curfman McInnes 
21659b94acceSBarry Smith    Input Parameter:
21669b94acceSBarry Smith .  snes - the SNES context
21679b94acceSBarry Smith 
21689b94acceSBarry Smith    Output Parameter:
21699b94acceSBarry Smith .  x - the solution
21709b94acceSBarry Smith 
217170e92668SMatthew Knepley    Level: intermediate
217236851e7fSLois Curfman McInnes 
21739b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
21749b94acceSBarry Smith 
217585385478SLisandro Dalcin .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
21769b94acceSBarry Smith @*/
217763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
21789b94acceSBarry Smith {
21793a40ed3dSBarry Smith   PetscFunctionBegin;
21804482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
21814482741eSBarry Smith   PetscValidPointer(x,2);
218285385478SLisandro Dalcin   *x = snes->vec_sol;
218370e92668SMatthew Knepley   PetscFunctionReturn(0);
218470e92668SMatthew Knepley }
218570e92668SMatthew Knepley 
218670e92668SMatthew Knepley #undef __FUNCT__
21874a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
218852baeb72SSatish Balay /*@
21899b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
21909b94acceSBarry Smith    stored.
21919b94acceSBarry Smith 
2192c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2193c7afd0dbSLois Curfman McInnes 
21949b94acceSBarry Smith    Input Parameter:
21959b94acceSBarry Smith .  snes - the SNES context
21969b94acceSBarry Smith 
21979b94acceSBarry Smith    Output Parameter:
21989b94acceSBarry Smith .  x - the solution update
21999b94acceSBarry Smith 
220036851e7fSLois Curfman McInnes    Level: advanced
220136851e7fSLois Curfman McInnes 
22029b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
22039b94acceSBarry Smith 
220485385478SLisandro Dalcin .seealso: SNESGetSolution(), SNESGetFunction()
22059b94acceSBarry Smith @*/
220663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
22079b94acceSBarry Smith {
22083a40ed3dSBarry Smith   PetscFunctionBegin;
22094482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
22104482741eSBarry Smith   PetscValidPointer(x,2);
221185385478SLisandro Dalcin   *x = snes->vec_sol_update;
22123a40ed3dSBarry Smith   PetscFunctionReturn(0);
22139b94acceSBarry Smith }
22149b94acceSBarry Smith 
22154a2ae208SSatish Balay #undef __FUNCT__
22164a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
22179b94acceSBarry Smith /*@C
22183638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
22199b94acceSBarry Smith 
2220c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2221c7afd0dbSLois Curfman McInnes 
22229b94acceSBarry Smith    Input Parameter:
22239b94acceSBarry Smith .  snes - the SNES context
22249b94acceSBarry Smith 
22259b94acceSBarry Smith    Output Parameter:
22267bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
222770e92668SMatthew Knepley .  func - the function (or PETSC_NULL)
222870e92668SMatthew Knepley -  ctx - the function context (or PETSC_NULL)
22299b94acceSBarry Smith 
223036851e7fSLois Curfman McInnes    Level: advanced
223136851e7fSLois Curfman McInnes 
2232a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
22339b94acceSBarry Smith 
22344b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
22359b94acceSBarry Smith @*/
223670e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
22379b94acceSBarry Smith {
22383a40ed3dSBarry Smith   PetscFunctionBegin;
22394482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
224085385478SLisandro Dalcin   if (r)    *r    = snes->vec_func;
2241e7788613SBarry Smith   if (func) *func = snes->ops->computefunction;
224270e92668SMatthew Knepley   if (ctx)  *ctx  = snes->funP;
22433a40ed3dSBarry Smith   PetscFunctionReturn(0);
22449b94acceSBarry Smith }
22459b94acceSBarry Smith 
22464a2ae208SSatish Balay #undef __FUNCT__
22474a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
22483c7409f5SSatish Balay /*@C
22493c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2250d850072dSLois Curfman McInnes    SNES options in the database.
22513c7409f5SSatish Balay 
2252fee21e36SBarry Smith    Collective on SNES
2253fee21e36SBarry Smith 
2254c7afd0dbSLois Curfman McInnes    Input Parameter:
2255c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2256c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2257c7afd0dbSLois Curfman McInnes 
2258d850072dSLois Curfman McInnes    Notes:
2259a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2260c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2261d850072dSLois Curfman McInnes 
226236851e7fSLois Curfman McInnes    Level: advanced
226336851e7fSLois Curfman McInnes 
22643c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2265a86d99e1SLois Curfman McInnes 
2266a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
22673c7409f5SSatish Balay @*/
226863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
22693c7409f5SSatish Balay {
2270dfbe8321SBarry Smith   PetscErrorCode ierr;
22713c7409f5SSatish Balay 
22723a40ed3dSBarry Smith   PetscFunctionBegin;
22734482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2274639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
227594b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
22763a40ed3dSBarry Smith   PetscFunctionReturn(0);
22773c7409f5SSatish Balay }
22783c7409f5SSatish Balay 
22794a2ae208SSatish Balay #undef __FUNCT__
22804a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
22813c7409f5SSatish Balay /*@C
2282f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2283d850072dSLois Curfman McInnes    SNES options in the database.
22843c7409f5SSatish Balay 
2285fee21e36SBarry Smith    Collective on SNES
2286fee21e36SBarry Smith 
2287c7afd0dbSLois Curfman McInnes    Input Parameters:
2288c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2289c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2290c7afd0dbSLois Curfman McInnes 
2291d850072dSLois Curfman McInnes    Notes:
2292a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2293c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2294d850072dSLois Curfman McInnes 
229536851e7fSLois Curfman McInnes    Level: advanced
229636851e7fSLois Curfman McInnes 
22973c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2298a86d99e1SLois Curfman McInnes 
2299a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
23003c7409f5SSatish Balay @*/
230163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
23023c7409f5SSatish Balay {
2303dfbe8321SBarry Smith   PetscErrorCode ierr;
23043c7409f5SSatish Balay 
23053a40ed3dSBarry Smith   PetscFunctionBegin;
23064482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2307639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
230894b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
23093a40ed3dSBarry Smith   PetscFunctionReturn(0);
23103c7409f5SSatish Balay }
23113c7409f5SSatish Balay 
23124a2ae208SSatish Balay #undef __FUNCT__
23134a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
23149ab63eb5SSatish Balay /*@C
23153c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
23163c7409f5SSatish Balay    SNES options in the database.
23173c7409f5SSatish Balay 
2318c7afd0dbSLois Curfman McInnes    Not Collective
2319c7afd0dbSLois Curfman McInnes 
23203c7409f5SSatish Balay    Input Parameter:
23213c7409f5SSatish Balay .  snes - the SNES context
23223c7409f5SSatish Balay 
23233c7409f5SSatish Balay    Output Parameter:
23243c7409f5SSatish Balay .  prefix - pointer to the prefix string used
23253c7409f5SSatish Balay 
23264ef407dbSRichard Tran Mills    Notes: On the fortran side, the user should pass in a string 'prefix' of
23279ab63eb5SSatish Balay    sufficient length to hold the prefix.
23289ab63eb5SSatish Balay 
232936851e7fSLois Curfman McInnes    Level: advanced
233036851e7fSLois Curfman McInnes 
23313c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2332a86d99e1SLois Curfman McInnes 
2333a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
23343c7409f5SSatish Balay @*/
2335e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
23363c7409f5SSatish Balay {
2337dfbe8321SBarry Smith   PetscErrorCode ierr;
23383c7409f5SSatish Balay 
23393a40ed3dSBarry Smith   PetscFunctionBegin;
23404482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2341639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
23423a40ed3dSBarry Smith   PetscFunctionReturn(0);
23433c7409f5SSatish Balay }
23443c7409f5SSatish Balay 
2345b2002411SLois Curfman McInnes 
23464a2ae208SSatish Balay #undef __FUNCT__
23474a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
23483cea93caSBarry Smith /*@C
23493cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
23503cea93caSBarry Smith 
23517f6c08e0SMatthew Knepley   Level: advanced
23523cea93caSBarry Smith @*/
235363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2354b2002411SLois Curfman McInnes {
2355e2d1d2b7SBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
2356dfbe8321SBarry Smith   PetscErrorCode ierr;
2357b2002411SLois Curfman McInnes 
2358b2002411SLois Curfman McInnes   PetscFunctionBegin;
2359b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2360c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2361b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2362b2002411SLois Curfman McInnes }
2363da9b6338SBarry Smith 
2364da9b6338SBarry Smith #undef __FUNCT__
2365da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
236663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2367da9b6338SBarry Smith {
2368dfbe8321SBarry Smith   PetscErrorCode ierr;
236977431f27SBarry Smith   PetscInt       N,i,j;
2370da9b6338SBarry Smith   Vec            u,uh,fh;
2371da9b6338SBarry Smith   PetscScalar    value;
2372da9b6338SBarry Smith   PetscReal      norm;
2373da9b6338SBarry Smith 
2374da9b6338SBarry Smith   PetscFunctionBegin;
2375da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2376da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2377da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2378da9b6338SBarry Smith 
2379da9b6338SBarry Smith   /* currently only works for sequential */
2380da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2381da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2382da9b6338SBarry Smith   for (i=0; i<N; i++) {
2383da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
238477431f27SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2385da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2386ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2387da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
23883ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2389da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
239077431f27SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2391da9b6338SBarry Smith       value = -value;
2392da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2393da9b6338SBarry Smith     }
2394da9b6338SBarry Smith   }
2395da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2396da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2397da9b6338SBarry Smith   PetscFunctionReturn(0);
2398da9b6338SBarry Smith }
239971f87433Sdalcinl 
240071f87433Sdalcinl #undef __FUNCT__
2401fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPSetUseEW"
240271f87433Sdalcinl /*@
2403fa9f3622SBarry Smith    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
240471f87433Sdalcinl    computing relative tolerance for linear solvers within an inexact
240571f87433Sdalcinl    Newton method.
240671f87433Sdalcinl 
240771f87433Sdalcinl    Collective on SNES
240871f87433Sdalcinl 
240971f87433Sdalcinl    Input Parameters:
241071f87433Sdalcinl +  snes - SNES context
241171f87433Sdalcinl -  flag - PETSC_TRUE or PETSC_FALSE
241271f87433Sdalcinl 
241371f87433Sdalcinl    Notes:
241471f87433Sdalcinl    Currently, the default is to use a constant relative tolerance for
241571f87433Sdalcinl    the inner linear solvers.  Alternatively, one can use the
241671f87433Sdalcinl    Eisenstat-Walker method, where the relative convergence tolerance
241771f87433Sdalcinl    is reset at each Newton iteration according progress of the nonlinear
241871f87433Sdalcinl    solver.
241971f87433Sdalcinl 
242071f87433Sdalcinl    Level: advanced
242171f87433Sdalcinl 
242271f87433Sdalcinl    Reference:
242371f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
242471f87433Sdalcinl    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
242571f87433Sdalcinl 
242671f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
242771f87433Sdalcinl 
2428fa9f3622SBarry Smith .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
242971f87433Sdalcinl @*/
2430fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag)
243171f87433Sdalcinl {
243271f87433Sdalcinl   PetscFunctionBegin;
243371f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
243471f87433Sdalcinl   snes->ksp_ewconv = flag;
243571f87433Sdalcinl   PetscFunctionReturn(0);
243671f87433Sdalcinl }
243771f87433Sdalcinl 
243871f87433Sdalcinl #undef __FUNCT__
2439fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPGetUseEW"
244071f87433Sdalcinl /*@
2441fa9f3622SBarry Smith    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
244271f87433Sdalcinl    for computing relative tolerance for linear solvers within an
244371f87433Sdalcinl    inexact Newton method.
244471f87433Sdalcinl 
244571f87433Sdalcinl    Not Collective
244671f87433Sdalcinl 
244771f87433Sdalcinl    Input Parameter:
244871f87433Sdalcinl .  snes - SNES context
244971f87433Sdalcinl 
245071f87433Sdalcinl    Output Parameter:
245171f87433Sdalcinl .  flag - PETSC_TRUE or PETSC_FALSE
245271f87433Sdalcinl 
245371f87433Sdalcinl    Notes:
245471f87433Sdalcinl    Currently, the default is to use a constant relative tolerance for
245571f87433Sdalcinl    the inner linear solvers.  Alternatively, one can use the
245671f87433Sdalcinl    Eisenstat-Walker method, where the relative convergence tolerance
245771f87433Sdalcinl    is reset at each Newton iteration according progress of the nonlinear
245871f87433Sdalcinl    solver.
245971f87433Sdalcinl 
246071f87433Sdalcinl    Level: advanced
246171f87433Sdalcinl 
246271f87433Sdalcinl    Reference:
246371f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
246471f87433Sdalcinl    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
246571f87433Sdalcinl 
246671f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
246771f87433Sdalcinl 
2468fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
246971f87433Sdalcinl @*/
2470fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag)
247171f87433Sdalcinl {
247271f87433Sdalcinl   PetscFunctionBegin;
247371f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
247471f87433Sdalcinl   PetscValidPointer(flag,2);
247571f87433Sdalcinl   *flag = snes->ksp_ewconv;
247671f87433Sdalcinl   PetscFunctionReturn(0);
247771f87433Sdalcinl }
247871f87433Sdalcinl 
247971f87433Sdalcinl #undef __FUNCT__
2480fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPSetParametersEW"
248171f87433Sdalcinl /*@
2482fa9f3622SBarry Smith    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
248371f87433Sdalcinl    convergence criteria for the linear solvers within an inexact
248471f87433Sdalcinl    Newton method.
248571f87433Sdalcinl 
248671f87433Sdalcinl    Collective on SNES
248771f87433Sdalcinl 
248871f87433Sdalcinl    Input Parameters:
248971f87433Sdalcinl +    snes - SNES context
249071f87433Sdalcinl .    version - version 1, 2 (default is 2) or 3
249171f87433Sdalcinl .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
249271f87433Sdalcinl .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
249371f87433Sdalcinl .    gamma - multiplicative factor for version 2 rtol computation
249471f87433Sdalcinl              (0 <= gamma2 <= 1)
249571f87433Sdalcinl .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
249671f87433Sdalcinl .    alpha2 - power for safeguard
249771f87433Sdalcinl -    threshold - threshold for imposing safeguard (0 < threshold < 1)
249871f87433Sdalcinl 
249971f87433Sdalcinl    Note:
250071f87433Sdalcinl    Version 3 was contributed by Luis Chacon, June 2006.
250171f87433Sdalcinl 
250271f87433Sdalcinl    Use PETSC_DEFAULT to retain the default for any of the parameters.
250371f87433Sdalcinl 
250471f87433Sdalcinl    Level: advanced
250571f87433Sdalcinl 
250671f87433Sdalcinl    Reference:
250771f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
250871f87433Sdalcinl    inexact Newton method", Utah State University Math. Stat. Dept. Res.
250971f87433Sdalcinl    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
251071f87433Sdalcinl 
251171f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
251271f87433Sdalcinl 
2513fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
251471f87433Sdalcinl @*/
2515fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
251671f87433Sdalcinl 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
251771f87433Sdalcinl {
2518fa9f3622SBarry Smith   SNESKSPEW *kctx;
251971f87433Sdalcinl   PetscFunctionBegin;
252071f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2521fa9f3622SBarry Smith   kctx = (SNESKSPEW*)snes->kspconvctx;
252271f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
252371f87433Sdalcinl 
252471f87433Sdalcinl   if (version != PETSC_DEFAULT)   kctx->version   = version;
252571f87433Sdalcinl   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
252671f87433Sdalcinl   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
252771f87433Sdalcinl   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
252871f87433Sdalcinl   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
252971f87433Sdalcinl   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
253071f87433Sdalcinl   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
253171f87433Sdalcinl 
253271f87433Sdalcinl   if (kctx->version < 1 || kctx->version > 3) {
253371f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
253471f87433Sdalcinl   }
253571f87433Sdalcinl   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
253671f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
253771f87433Sdalcinl   }
253871f87433Sdalcinl   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
253971f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
254071f87433Sdalcinl   }
254171f87433Sdalcinl   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
254271f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
254371f87433Sdalcinl   }
254471f87433Sdalcinl   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
254571f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
254671f87433Sdalcinl   }
254771f87433Sdalcinl   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
254871f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
254971f87433Sdalcinl   }
255071f87433Sdalcinl   PetscFunctionReturn(0);
255171f87433Sdalcinl }
255271f87433Sdalcinl 
255371f87433Sdalcinl #undef __FUNCT__
2554fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPGetParametersEW"
255571f87433Sdalcinl /*@
2556fa9f3622SBarry Smith    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
255771f87433Sdalcinl    convergence criteria for the linear solvers within an inexact
255871f87433Sdalcinl    Newton method.
255971f87433Sdalcinl 
256071f87433Sdalcinl    Not Collective
256171f87433Sdalcinl 
256271f87433Sdalcinl    Input Parameters:
256371f87433Sdalcinl      snes - SNES context
256471f87433Sdalcinl 
256571f87433Sdalcinl    Output Parameters:
256671f87433Sdalcinl +    version - version 1, 2 (default is 2) or 3
256771f87433Sdalcinl .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
256871f87433Sdalcinl .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
256971f87433Sdalcinl .    gamma - multiplicative factor for version 2 rtol computation
257071f87433Sdalcinl              (0 <= gamma2 <= 1)
257171f87433Sdalcinl .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
257271f87433Sdalcinl .    alpha2 - power for safeguard
257371f87433Sdalcinl -    threshold - threshold for imposing safeguard (0 < threshold < 1)
257471f87433Sdalcinl 
257571f87433Sdalcinl    Level: advanced
257671f87433Sdalcinl 
257771f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
257871f87433Sdalcinl 
2579fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
258071f87433Sdalcinl @*/
2581fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
258271f87433Sdalcinl 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
258371f87433Sdalcinl {
2584fa9f3622SBarry Smith   SNESKSPEW *kctx;
258571f87433Sdalcinl   PetscFunctionBegin;
258671f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2587fa9f3622SBarry Smith   kctx = (SNESKSPEW*)snes->kspconvctx;
258871f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
258971f87433Sdalcinl   if(version)   *version   = kctx->version;
259071f87433Sdalcinl   if(rtol_0)    *rtol_0    = kctx->rtol_0;
259171f87433Sdalcinl   if(rtol_max)  *rtol_max  = kctx->rtol_max;
259271f87433Sdalcinl   if(gamma)     *gamma     = kctx->gamma;
259371f87433Sdalcinl   if(alpha)     *alpha     = kctx->alpha;
259471f87433Sdalcinl   if(alpha2)    *alpha2    = kctx->alpha2;
259571f87433Sdalcinl   if(threshold) *threshold = kctx->threshold;
259671f87433Sdalcinl   PetscFunctionReturn(0);
259771f87433Sdalcinl }
259871f87433Sdalcinl 
259971f87433Sdalcinl #undef __FUNCT__
2600fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPEW_PreSolve"
2601fa9f3622SBarry Smith static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
260271f87433Sdalcinl {
260371f87433Sdalcinl   PetscErrorCode ierr;
2604fa9f3622SBarry Smith   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
260571f87433Sdalcinl   PetscReal      rtol=PETSC_DEFAULT,stol;
260671f87433Sdalcinl 
260771f87433Sdalcinl   PetscFunctionBegin;
260871f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
260971f87433Sdalcinl   if (!snes->iter) { /* first time in, so use the original user rtol */
261071f87433Sdalcinl     rtol = kctx->rtol_0;
261171f87433Sdalcinl   } else {
261271f87433Sdalcinl     if (kctx->version == 1) {
261371f87433Sdalcinl       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
261471f87433Sdalcinl       if (rtol < 0.0) rtol = -rtol;
261571f87433Sdalcinl       stol = pow(kctx->rtol_last,kctx->alpha2);
261671f87433Sdalcinl       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
261771f87433Sdalcinl     } else if (kctx->version == 2) {
261871f87433Sdalcinl       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
261971f87433Sdalcinl       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
262071f87433Sdalcinl       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
262171f87433Sdalcinl     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
262271f87433Sdalcinl       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
262371f87433Sdalcinl       /* safeguard: avoid sharp decrease of rtol */
262471f87433Sdalcinl       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
262571f87433Sdalcinl       stol = PetscMax(rtol,stol);
262671f87433Sdalcinl       rtol = PetscMin(kctx->rtol_0,stol);
262771f87433Sdalcinl       /* safeguard: avoid oversolving */
262871f87433Sdalcinl       stol = kctx->gamma*(snes->ttol)/snes->norm;
262971f87433Sdalcinl       stol = PetscMax(rtol,stol);
263071f87433Sdalcinl       rtol = PetscMin(kctx->rtol_0,stol);
263171f87433Sdalcinl     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
263271f87433Sdalcinl   }
263371f87433Sdalcinl   /* safeguard: avoid rtol greater than one */
263471f87433Sdalcinl   rtol = PetscMin(rtol,kctx->rtol_max);
263571f87433Sdalcinl   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
263671f87433Sdalcinl   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
263771f87433Sdalcinl   PetscFunctionReturn(0);
263871f87433Sdalcinl }
263971f87433Sdalcinl 
264071f87433Sdalcinl #undef __FUNCT__
2641fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPEW_PostSolve"
2642fa9f3622SBarry Smith static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
264371f87433Sdalcinl {
264471f87433Sdalcinl   PetscErrorCode ierr;
2645fa9f3622SBarry Smith   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
264671f87433Sdalcinl   PCSide         pcside;
264771f87433Sdalcinl   Vec            lres;
264871f87433Sdalcinl 
264971f87433Sdalcinl   PetscFunctionBegin;
265071f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
265171f87433Sdalcinl   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
265271f87433Sdalcinl   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
265371f87433Sdalcinl   if (kctx->version == 1) {
265471f87433Sdalcinl     ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr);
265571f87433Sdalcinl     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
265671f87433Sdalcinl       /* KSP residual is true linear residual */
265771f87433Sdalcinl       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
265871f87433Sdalcinl     } else {
265971f87433Sdalcinl       /* KSP residual is preconditioned residual */
266071f87433Sdalcinl       /* compute true linear residual norm */
266171f87433Sdalcinl       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
266271f87433Sdalcinl       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
266371f87433Sdalcinl       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
266471f87433Sdalcinl       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
266571f87433Sdalcinl       ierr = VecDestroy(lres);CHKERRQ(ierr);
266671f87433Sdalcinl     }
266771f87433Sdalcinl   }
266871f87433Sdalcinl   PetscFunctionReturn(0);
266971f87433Sdalcinl }
267071f87433Sdalcinl 
267171f87433Sdalcinl #undef __FUNCT__
267271f87433Sdalcinl #define __FUNCT__ "SNES_KSPSolve"
267371f87433Sdalcinl PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
267471f87433Sdalcinl {
267571f87433Sdalcinl   PetscErrorCode ierr;
267671f87433Sdalcinl 
267771f87433Sdalcinl   PetscFunctionBegin;
2678fa9f3622SBarry Smith   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
267971f87433Sdalcinl   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
2680fa9f3622SBarry Smith   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
268171f87433Sdalcinl   PetscFunctionReturn(0);
268271f87433Sdalcinl }
2683