xref: /petsc/src/snes/interface/snes.c (revision 58ebbce7ab732c9f5ad3d7b5906f4363e35c43a5)
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 */
963dd3a1aSKris Buschelman PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE = 0;
106849ba73SBarry Smith PetscEvent  SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0;
11a09944afSBarry Smith 
12a09944afSBarry Smith #undef __FUNCT__
134a2ae208SSatish Balay #define __FUNCT__ "SNESView"
147e2c5f70SBarry Smith /*@C
159b94acceSBarry Smith    SNESView - Prints the SNES data structure.
169b94acceSBarry Smith 
174c49b128SBarry Smith    Collective on SNES
18fee21e36SBarry Smith 
19c7afd0dbSLois Curfman McInnes    Input Parameters:
20c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
21c7afd0dbSLois Curfman McInnes -  viewer - visualization context
22c7afd0dbSLois Curfman McInnes 
239b94acceSBarry Smith    Options Database Key:
24c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
259b94acceSBarry Smith 
269b94acceSBarry Smith    Notes:
279b94acceSBarry Smith    The available visualization contexts include
28b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
29b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
30c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
31c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
32c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
339b94acceSBarry Smith 
343e081fefSLois Curfman McInnes    The user can open an alternative visualization context with
35b0a32e0cSBarry Smith    PetscViewerASCIIOpen() - output to a specified file.
369b94acceSBarry Smith 
3736851e7fSLois Curfman McInnes    Level: beginner
3836851e7fSLois Curfman McInnes 
399b94acceSBarry Smith .keywords: SNES, view
409b94acceSBarry Smith 
41b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
429b94acceSBarry Smith @*/
4363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer)
449b94acceSBarry Smith {
45fa9f3622SBarry Smith   SNESKSPEW           *kctx;
46dfbe8321SBarry Smith   PetscErrorCode      ierr;
4794b7f48cSBarry Smith   KSP                 ksp;
48e060cb09SBarry Smith   SNESType            type;
4932077d6dSBarry Smith   PetscTruth          iascii,isstring;
509b94acceSBarry Smith 
513a40ed3dSBarry Smith   PetscFunctionBegin;
524482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
533050cee2SBarry Smith   if (!viewer) {
547adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr);
553050cee2SBarry Smith   }
564482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
57c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,viewer,2);
5874679c65SBarry Smith 
5932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
60b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
6132077d6dSBarry Smith   if (iascii) {
627adad957SLisandro Dalcin     if (((PetscObject)snes)->prefix) {
637adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr);
643a7fca6bSBarry Smith     } else {
65b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
663a7fca6bSBarry Smith     }
67454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
68454a90a3SBarry Smith     if (type) {
69b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
70184914b5SBarry Smith     } else {
71b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
72184914b5SBarry Smith     }
73e7788613SBarry Smith     if (snes->ops->view) {
74b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
75e7788613SBarry Smith       ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr);
76b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
770ef38995SBarry Smith     }
7877431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
79a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
8070441072SBarry Smith                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
8177431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
8277431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
839b94acceSBarry Smith     if (snes->ksp_ewconv) {
84fa9f3622SBarry Smith       kctx = (SNESKSPEW *)snes->kspconvctx;
859b94acceSBarry Smith       if (kctx) {
8677431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
87a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
88a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
899b94acceSBarry Smith       }
909b94acceSBarry Smith     }
910f5bd95cSBarry Smith   } else if (isstring) {
92454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
93b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
9419bcc07fSBarry Smith   }
9594b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
96b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
9794b7f48cSBarry Smith   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
98b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
993a40ed3dSBarry Smith   PetscFunctionReturn(0);
1009b94acceSBarry Smith }
1019b94acceSBarry Smith 
10276b2cf59SMatthew Knepley /*
10376b2cf59SMatthew Knepley   We retain a list of functions that also take SNES command
10476b2cf59SMatthew Knepley   line options. These are called at the end SNESSetFromOptions()
10576b2cf59SMatthew Knepley */
10676b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5
107a7cc72afSBarry Smith static PetscInt numberofsetfromoptions;
1086849ba73SBarry Smith static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
10976b2cf59SMatthew Knepley 
110e74ef692SMatthew Knepley #undef __FUNCT__
111e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker"
112ac226902SBarry Smith /*@C
11376b2cf59SMatthew Knepley   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
11476b2cf59SMatthew Knepley 
11576b2cf59SMatthew Knepley   Not Collective
11676b2cf59SMatthew Knepley 
11776b2cf59SMatthew Knepley   Input Parameter:
11876b2cf59SMatthew Knepley . snescheck - function that checks for options
11976b2cf59SMatthew Knepley 
12076b2cf59SMatthew Knepley   Level: developer
12176b2cf59SMatthew Knepley 
12276b2cf59SMatthew Knepley .seealso: SNESSetFromOptions()
12376b2cf59SMatthew Knepley @*/
12463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
12576b2cf59SMatthew Knepley {
12676b2cf59SMatthew Knepley   PetscFunctionBegin;
12776b2cf59SMatthew Knepley   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
12877431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
12976b2cf59SMatthew Knepley   }
13076b2cf59SMatthew Knepley   othersetfromoptions[numberofsetfromoptions++] = snescheck;
13176b2cf59SMatthew Knepley   PetscFunctionReturn(0);
13276b2cf59SMatthew Knepley }
13376b2cf59SMatthew Knepley 
1344a2ae208SSatish Balay #undef __FUNCT__
1354a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions"
1369b94acceSBarry Smith /*@
13794b7f48cSBarry Smith    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
1389b94acceSBarry Smith 
139c7afd0dbSLois Curfman McInnes    Collective on SNES
140c7afd0dbSLois Curfman McInnes 
1419b94acceSBarry Smith    Input Parameter:
1429b94acceSBarry Smith .  snes - the SNES context
1439b94acceSBarry Smith 
14436851e7fSLois Curfman McInnes    Options Database Keys:
1456831982aSBarry Smith +  -snes_type <type> - ls, tr, umls, umtr, test
14682738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
14782738288SBarry Smith                 of the change in the solution between steps
14870441072SBarry Smith .  -snes_atol <abstol> - absolute tolerance of residual norm
149b39c3a46SLois Curfman McInnes .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
150b39c3a46SLois Curfman McInnes .  -snes_max_it <max_it> - maximum number of iterations
151b39c3a46SLois Curfman McInnes .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
15250ffb88aSMatthew Knepley .  -snes_max_fail <max_fail> - maximum number of failures
153b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
1542492ecdbSBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear
15582738288SBarry Smith                                solver; hence iterations will continue until max_it
1561fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
15782738288SBarry Smith                                of convergence test
158e8105e01SRichard Katz .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
159e8105e01SRichard Katz                                        filename given prints to stdout
160a6570f20SBarry Smith .  -snes_monitor_solution - plots solution at each iteration
161a6570f20SBarry Smith .  -snes_monitor_residual - plots residual (not its norm) at each iteration
162a6570f20SBarry Smith .  -snes_monitor_solution_update - plots update to solution at each iteration
163a6570f20SBarry Smith .  -snes_monitor_draw - plots residual norm at each iteration
164e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
1655968eb51SBarry Smith .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
166fee2055bSBarry Smith -  -snes_converged_reason - print the reason for convergence/divergence after each solve
16782738288SBarry Smith 
16882738288SBarry Smith     Options Database for Eisenstat-Walker method:
169fa9f3622SBarry Smith +  -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence
1704b27c08aSLois Curfman McInnes .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
17136851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
17236851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
17336851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
17436851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
17536851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
17636851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
17782738288SBarry Smith 
17811ca99fdSLois Curfman McInnes    Notes:
17911ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
18011ca99fdSLois Curfman McInnes    the users manual.
18183e2fdc7SBarry Smith 
18236851e7fSLois Curfman McInnes    Level: beginner
18336851e7fSLois Curfman McInnes 
1849b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1859b94acceSBarry Smith 
18669ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
1879b94acceSBarry Smith @*/
18863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes)
1899b94acceSBarry Smith {
190f1af5d2fSBarry Smith   PetscTruth              flg;
19185385478SLisandro Dalcin   PetscInt                i,indx;
19285385478SLisandro Dalcin   const char              *deft = SNESLS;
19385385478SLisandro Dalcin   const char              *convtests[] = {"default","skip"};
19485385478SLisandro Dalcin   SNESKSPEW               *kctx = NULL;
195e8105e01SRichard Katz   char                    type[256], monfilename[PETSC_MAX_PATH_LEN];
19623d894e5SBarry Smith   PetscViewerASCIIMonitor monviewer;
19785385478SLisandro Dalcin   PetscErrorCode          ierr;
1989b94acceSBarry Smith 
1993a40ed3dSBarry Smith   PetscFunctionBegin;
2004482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
201ca161407SBarry Smith 
2027adad957SLisandro Dalcin   ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
203186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
2047adad957SLisandro Dalcin     if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; }
205b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
206d64ed03dSBarry Smith     if (flg) {
207186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
2087adad957SLisandro Dalcin     } else if (!((PetscObject)snes)->type_name) {
209186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
210d64ed03dSBarry Smith     }
211909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
21293c39befSBarry Smith 
21387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
21470441072SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
215186905e3SBarry Smith 
21687828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
217b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
218b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
21950ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
22085385478SLisandro Dalcin 
22185385478SLisandro Dalcin     ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr);
22285385478SLisandro Dalcin     if (flg) {
22385385478SLisandro Dalcin       switch (indx) {
22485385478SLisandro Dalcin       case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL);CHKERRQ(ierr); break;
22585385478SLisandro Dalcin       case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL);CHKERRQ(ierr);    break;
22685385478SLisandro Dalcin       }
22785385478SLisandro Dalcin     }
22885385478SLisandro Dalcin 
2295968eb51SBarry Smith     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
2305968eb51SBarry Smith     if (flg) {
2315968eb51SBarry Smith       snes->printreason = PETSC_TRUE;
2325968eb51SBarry Smith     }
233186905e3SBarry Smith 
23485385478SLisandro Dalcin     kctx = (SNESKSPEW *)snes->kspconvctx;
23585385478SLisandro Dalcin 
236fa9f3622SBarry Smith     ierr = PetscOptionsTruth("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr);
237186905e3SBarry Smith 
238fa9f3622SBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
239fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
240fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
241fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
242fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
243fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
244fa9f3622SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
245186905e3SBarry Smith 
246a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",&flg);CHKERRQ(ierr);
247a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);}
248eabae89aSBarry Smith 
249a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
250e8105e01SRichard Katz     if (flg) {
2517adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
25223d894e5SBarry Smith       ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
253e8105e01SRichard Katz     }
254eabae89aSBarry Smith 
255a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
256eabae89aSBarry Smith     if (flg) {
2577adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
258f1bef1bcSMatthew Knepley       ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr);
259e8105e01SRichard Katz     }
260eabae89aSBarry Smith 
261a6570f20SBarry Smith     ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
262eabae89aSBarry Smith     if (flg) {
2637adad957SLisandro Dalcin       ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,0,&monviewer);CHKERRQ(ierr);
26423d894e5SBarry Smith       ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
265eabae89aSBarry Smith     }
266eabae89aSBarry Smith 
267a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",&flg);CHKERRQ(ierr);
268a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);}
269a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",&flg);CHKERRQ(ierr);
270a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);}
271a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",&flg);CHKERRQ(ierr);
272a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);}
273a6570f20SBarry Smith     ierr = PetscOptionsName("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",&flg);CHKERRQ(ierr);
274a6570f20SBarry Smith     if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
275e24b481bSBarry Smith 
276b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2774b27c08aSLois Curfman McInnes     if (flg) {
278186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
279ae15b995SBarry Smith       ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr);
2809b94acceSBarry Smith     }
281639f9d9dSBarry Smith 
28276b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
28376b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
28476b2cf59SMatthew Knepley     }
28576b2cf59SMatthew Knepley 
286e7788613SBarry Smith     if (snes->ops->setfromoptions) {
287e7788613SBarry Smith       ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr);
288639f9d9dSBarry Smith     }
289b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2904bbc92c1SBarry Smith 
29185385478SLisandro Dalcin   ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr);
29293993e2dSLois Curfman McInnes 
2933a40ed3dSBarry Smith   PetscFunctionReturn(0);
2949b94acceSBarry Smith }
2959b94acceSBarry Smith 
296a847f771SSatish Balay 
2974a2ae208SSatish Balay #undef __FUNCT__
2984a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
2999b94acceSBarry Smith /*@
3009b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3019b94acceSBarry Smith    the nonlinear solvers.
3029b94acceSBarry Smith 
303fee21e36SBarry Smith    Collective on SNES
304fee21e36SBarry Smith 
305c7afd0dbSLois Curfman McInnes    Input Parameters:
306c7afd0dbSLois Curfman McInnes +  snes - the SNES context
307c7afd0dbSLois Curfman McInnes -  usrP - optional user context
308c7afd0dbSLois Curfman McInnes 
30936851e7fSLois Curfman McInnes    Level: intermediate
31036851e7fSLois Curfman McInnes 
3119b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3129b94acceSBarry Smith 
3139b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3149b94acceSBarry Smith @*/
31563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP)
3169b94acceSBarry Smith {
3173a40ed3dSBarry Smith   PetscFunctionBegin;
3184482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3199b94acceSBarry Smith   snes->user		= usrP;
3203a40ed3dSBarry Smith   PetscFunctionReturn(0);
3219b94acceSBarry Smith }
32274679c65SBarry Smith 
3234a2ae208SSatish Balay #undef __FUNCT__
3244a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3259b94acceSBarry Smith /*@C
3269b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3279b94acceSBarry Smith    nonlinear solvers.
3289b94acceSBarry Smith 
329c7afd0dbSLois Curfman McInnes    Not Collective
330c7afd0dbSLois Curfman McInnes 
3319b94acceSBarry Smith    Input Parameter:
3329b94acceSBarry Smith .  snes - SNES context
3339b94acceSBarry Smith 
3349b94acceSBarry Smith    Output Parameter:
3359b94acceSBarry Smith .  usrP - user context
3369b94acceSBarry Smith 
33736851e7fSLois Curfman McInnes    Level: intermediate
33836851e7fSLois Curfman McInnes 
3399b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3409b94acceSBarry Smith 
3419b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3429b94acceSBarry Smith @*/
34363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP)
3449b94acceSBarry Smith {
3453a40ed3dSBarry Smith   PetscFunctionBegin;
3464482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3479b94acceSBarry Smith   *usrP = snes->user;
3483a40ed3dSBarry Smith   PetscFunctionReturn(0);
3499b94acceSBarry Smith }
35074679c65SBarry Smith 
3514a2ae208SSatish Balay #undef __FUNCT__
3524a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3539b94acceSBarry Smith /*@
354c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
355c8228a4eSBarry Smith    at this time.
3569b94acceSBarry Smith 
357c7afd0dbSLois Curfman McInnes    Not Collective
358c7afd0dbSLois Curfman McInnes 
3599b94acceSBarry Smith    Input Parameter:
3609b94acceSBarry Smith .  snes - SNES context
3619b94acceSBarry Smith 
3629b94acceSBarry Smith    Output Parameter:
3639b94acceSBarry Smith .  iter - iteration number
3649b94acceSBarry Smith 
365c8228a4eSBarry Smith    Notes:
366c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
367c8228a4eSBarry Smith 
368c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
36908405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
37008405cd6SLois Curfman McInnes .vb
37108405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
37208405cd6SLois Curfman McInnes       if (!(it % 2)) {
37308405cd6SLois Curfman McInnes         [compute Jacobian here]
37408405cd6SLois Curfman McInnes       }
37508405cd6SLois Curfman McInnes .ve
376c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
37708405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
378c8228a4eSBarry Smith 
37936851e7fSLois Curfman McInnes    Level: intermediate
38036851e7fSLois Curfman McInnes 
3812b668275SBarry Smith .keywords: SNES, nonlinear, get, iteration, number,
3822b668275SBarry Smith 
383b850b91aSLisandro Dalcin .seealso:   SNESGetFunctionNorm(), SNESGetLinearSolveIterations()
3849b94acceSBarry Smith @*/
38563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter)
3869b94acceSBarry Smith {
3873a40ed3dSBarry Smith   PetscFunctionBegin;
3884482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3894482741eSBarry Smith   PetscValidIntPointer(iter,2);
3909b94acceSBarry Smith   *iter = snes->iter;
3913a40ed3dSBarry Smith   PetscFunctionReturn(0);
3929b94acceSBarry Smith }
39374679c65SBarry Smith 
3944a2ae208SSatish Balay #undef __FUNCT__
3954a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
3969b94acceSBarry Smith /*@
3979b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3989b94acceSBarry Smith    with SNESSSetFunction().
3999b94acceSBarry Smith 
400c7afd0dbSLois Curfman McInnes    Collective on SNES
401c7afd0dbSLois Curfman McInnes 
4029b94acceSBarry Smith    Input Parameter:
4039b94acceSBarry Smith .  snes - SNES context
4049b94acceSBarry Smith 
4059b94acceSBarry Smith    Output Parameter:
4069b94acceSBarry Smith .  fnorm - 2-norm of function
4079b94acceSBarry Smith 
40836851e7fSLois Curfman McInnes    Level: intermediate
40936851e7fSLois Curfman McInnes 
4109b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
411a86d99e1SLois Curfman McInnes 
412b850b91aSLisandro Dalcin .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations()
4139b94acceSBarry Smith @*/
41471f87433Sdalcinl PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscReal *fnorm)
4159b94acceSBarry Smith {
4163a40ed3dSBarry Smith   PetscFunctionBegin;
4174482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4184482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
4199b94acceSBarry Smith   *fnorm = snes->norm;
4203a40ed3dSBarry Smith   PetscFunctionReturn(0);
4219b94acceSBarry Smith }
42274679c65SBarry Smith 
4234a2ae208SSatish Balay #undef __FUNCT__
424b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetNonlinearStepFailures"
4259b94acceSBarry Smith /*@
426b850b91aSLisandro Dalcin    SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps
4279b94acceSBarry Smith    attempted by the nonlinear solver.
4289b94acceSBarry Smith 
429c7afd0dbSLois Curfman McInnes    Not Collective
430c7afd0dbSLois Curfman McInnes 
4319b94acceSBarry Smith    Input Parameter:
4329b94acceSBarry Smith .  snes - SNES context
4339b94acceSBarry Smith 
4349b94acceSBarry Smith    Output Parameter:
4359b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4369b94acceSBarry Smith 
437c96a6f78SLois Curfman McInnes    Notes:
438c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
439c96a6f78SLois Curfman McInnes 
44036851e7fSLois Curfman McInnes    Level: intermediate
44136851e7fSLois Curfman McInnes 
4429b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
443*58ebbce7SBarry Smith 
444*58ebbce7SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolverIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
445*58ebbce7SBarry Smith           SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures()
4469b94acceSBarry Smith @*/
447b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails)
4489b94acceSBarry Smith {
4493a40ed3dSBarry Smith   PetscFunctionBegin;
4504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4514482741eSBarry Smith   PetscValidIntPointer(nfails,2);
45250ffb88aSMatthew Knepley   *nfails = snes->numFailures;
45350ffb88aSMatthew Knepley   PetscFunctionReturn(0);
45450ffb88aSMatthew Knepley }
45550ffb88aSMatthew Knepley 
45650ffb88aSMatthew Knepley #undef __FUNCT__
457b850b91aSLisandro Dalcin #define __FUNCT__ "SNESSetMaxNonlinearStepFailures"
45850ffb88aSMatthew Knepley /*@
459b850b91aSLisandro Dalcin    SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps
46050ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
46150ffb88aSMatthew Knepley 
46250ffb88aSMatthew Knepley    Not Collective
46350ffb88aSMatthew Knepley 
46450ffb88aSMatthew Knepley    Input Parameters:
46550ffb88aSMatthew Knepley +  snes     - SNES context
46650ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
46750ffb88aSMatthew Knepley 
46850ffb88aSMatthew Knepley    Level: intermediate
46950ffb88aSMatthew Knepley 
47050ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
471*58ebbce7SBarry Smith 
472*58ebbce7SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolverIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
473*58ebbce7SBarry Smith           SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
47450ffb88aSMatthew Knepley @*/
475b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails)
47650ffb88aSMatthew Knepley {
47750ffb88aSMatthew Knepley   PetscFunctionBegin;
4784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
47950ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
48050ffb88aSMatthew Knepley   PetscFunctionReturn(0);
48150ffb88aSMatthew Knepley }
48250ffb88aSMatthew Knepley 
48350ffb88aSMatthew Knepley #undef __FUNCT__
484b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetMaxNonlinearStepFailures"
48550ffb88aSMatthew Knepley /*@
486b850b91aSLisandro Dalcin    SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps
48750ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
48850ffb88aSMatthew Knepley 
48950ffb88aSMatthew Knepley    Not Collective
49050ffb88aSMatthew Knepley 
49150ffb88aSMatthew Knepley    Input Parameter:
49250ffb88aSMatthew Knepley .  snes     - SNES context
49350ffb88aSMatthew Knepley 
49450ffb88aSMatthew Knepley    Output Parameter:
49550ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
49650ffb88aSMatthew Knepley 
49750ffb88aSMatthew Knepley    Level: intermediate
49850ffb88aSMatthew Knepley 
49950ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
500*58ebbce7SBarry Smith 
501*58ebbce7SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolverIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(),
502*58ebbce7SBarry Smith           SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures()
503*58ebbce7SBarry Smith 
50450ffb88aSMatthew Knepley @*/
505b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails)
50650ffb88aSMatthew Knepley {
50750ffb88aSMatthew Knepley   PetscFunctionBegin;
5084482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5094482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
51050ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
5113a40ed3dSBarry Smith   PetscFunctionReturn(0);
5129b94acceSBarry Smith }
513a847f771SSatish Balay 
5144a2ae208SSatish Balay #undef __FUNCT__
5152541af92SBarry Smith #define __FUNCT__ "SNESGetNumberFunctionEvals"
5162541af92SBarry Smith /*@
5172541af92SBarry Smith    SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations
5182541af92SBarry Smith      done by SNES.
5192541af92SBarry Smith 
5202541af92SBarry Smith    Not Collective
5212541af92SBarry Smith 
5222541af92SBarry Smith    Input Parameter:
5232541af92SBarry Smith .  snes     - SNES context
5242541af92SBarry Smith 
5252541af92SBarry Smith    Output Parameter:
5262541af92SBarry Smith .  nfuncs - number of evaluations
5272541af92SBarry Smith 
5282541af92SBarry Smith    Level: intermediate
5292541af92SBarry Smith 
5302541af92SBarry Smith .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
531*58ebbce7SBarry Smith 
532*58ebbce7SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolverIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures()
5332541af92SBarry Smith @*/
5342541af92SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs)
5352541af92SBarry Smith {
5362541af92SBarry Smith   PetscFunctionBegin;
5372541af92SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5382541af92SBarry Smith   PetscValidIntPointer(nfuncs,2);
5392541af92SBarry Smith   *nfuncs = snes->nfuncs;
5402541af92SBarry Smith   PetscFunctionReturn(0);
5412541af92SBarry Smith }
5422541af92SBarry Smith 
5432541af92SBarry Smith #undef __FUNCT__
5443d4c4710SBarry Smith #define __FUNCT__ "SNESGetLinearSolveFailures"
5453d4c4710SBarry Smith /*@
5463d4c4710SBarry Smith    SNESGetLinearSolveFailures - Gets the number of failed (non-converged)
5473d4c4710SBarry Smith    linear solvers.
5483d4c4710SBarry Smith 
5493d4c4710SBarry Smith    Not Collective
5503d4c4710SBarry Smith 
5513d4c4710SBarry Smith    Input Parameter:
5523d4c4710SBarry Smith .  snes - SNES context
5533d4c4710SBarry Smith 
5543d4c4710SBarry Smith    Output Parameter:
5553d4c4710SBarry Smith .  nfails - number of failed solves
5563d4c4710SBarry Smith 
5573d4c4710SBarry Smith    Notes:
5583d4c4710SBarry Smith    This counter is reset to zero for each successive call to SNESSolve().
5593d4c4710SBarry Smith 
5603d4c4710SBarry Smith    Level: intermediate
5613d4c4710SBarry Smith 
5623d4c4710SBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
563*58ebbce7SBarry Smith 
564*58ebbce7SBarry Smith .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolverIterations(), SNESSetMaxLinearSolveFailures()
5653d4c4710SBarry Smith @*/
5663d4c4710SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails)
5673d4c4710SBarry Smith {
5683d4c4710SBarry Smith   PetscFunctionBegin;
5693d4c4710SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5703d4c4710SBarry Smith   PetscValidIntPointer(nfails,2);
5713d4c4710SBarry Smith   *nfails = snes->numLinearSolveFailures;
5723d4c4710SBarry Smith   PetscFunctionReturn(0);
5733d4c4710SBarry Smith }
5743d4c4710SBarry Smith 
5753d4c4710SBarry Smith #undef __FUNCT__
5763d4c4710SBarry Smith #define __FUNCT__ "SNESSetMaxLinearSolveFailures"
5773d4c4710SBarry Smith /*@
5783d4c4710SBarry Smith    SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts
5793d4c4710SBarry Smith    allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE
5803d4c4710SBarry Smith 
5813d4c4710SBarry Smith    Collective on SNES
5823d4c4710SBarry Smith 
5833d4c4710SBarry Smith    Input Parameters:
5843d4c4710SBarry Smith +  snes     - SNES context
5853d4c4710SBarry Smith -  maxFails - maximum allowed linear solve failures
5863d4c4710SBarry Smith 
5873d4c4710SBarry Smith    Level: intermediate
5883d4c4710SBarry Smith 
5893d4c4710SBarry Smith    Notes: By default this is 1; that is SNES returns on the first failed linear solve
5903d4c4710SBarry Smith 
5913d4c4710SBarry Smith .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
5923d4c4710SBarry Smith 
593*58ebbce7SBarry Smith .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations()
5943d4c4710SBarry Smith @*/
5953d4c4710SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails)
5963d4c4710SBarry Smith {
5973d4c4710SBarry Smith   PetscFunctionBegin;
5983d4c4710SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5993d4c4710SBarry Smith   snes->maxLinearSolveFailures = maxFails;
6003d4c4710SBarry Smith   PetscFunctionReturn(0);
6013d4c4710SBarry Smith }
6023d4c4710SBarry Smith 
6033d4c4710SBarry Smith #undef __FUNCT__
6043d4c4710SBarry Smith #define __FUNCT__ "SNESGetMaxLinearSolveFailures"
6053d4c4710SBarry Smith /*@
6063d4c4710SBarry Smith    SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that
6073d4c4710SBarry Smith      are allowed before SNES terminates
6083d4c4710SBarry Smith 
6093d4c4710SBarry Smith    Not Collective
6103d4c4710SBarry Smith 
6113d4c4710SBarry Smith    Input Parameter:
6123d4c4710SBarry Smith .  snes     - SNES context
6133d4c4710SBarry Smith 
6143d4c4710SBarry Smith    Output Parameter:
6153d4c4710SBarry Smith .  maxFails - maximum of unsuccessful solves allowed
6163d4c4710SBarry Smith 
6173d4c4710SBarry Smith    Level: intermediate
6183d4c4710SBarry Smith 
6193d4c4710SBarry Smith    Notes: By default this is 1; that is SNES returns on the first failed linear solve
6203d4c4710SBarry Smith 
6213d4c4710SBarry Smith .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
6223d4c4710SBarry Smith 
623*58ebbce7SBarry Smith .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolverIterations(), SNESSetMaxLinearSolveFailures(),
6243d4c4710SBarry Smith @*/
6253d4c4710SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails)
6263d4c4710SBarry Smith {
6273d4c4710SBarry Smith   PetscFunctionBegin;
6283d4c4710SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6293d4c4710SBarry Smith   PetscValidIntPointer(maxFails,2);
6303d4c4710SBarry Smith   *maxFails = snes->maxLinearSolveFailures;
6313d4c4710SBarry Smith   PetscFunctionReturn(0);
6323d4c4710SBarry Smith }
6333d4c4710SBarry Smith 
6343d4c4710SBarry Smith #undef __FUNCT__
635b850b91aSLisandro Dalcin #define __FUNCT__ "SNESGetLinearSolveIterations"
636c96a6f78SLois Curfman McInnes /*@
637b850b91aSLisandro Dalcin    SNESGetLinearSolveIterations - Gets the total number of linear iterations
638c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
639c96a6f78SLois Curfman McInnes 
640c7afd0dbSLois Curfman McInnes    Not Collective
641c7afd0dbSLois Curfman McInnes 
642c96a6f78SLois Curfman McInnes    Input Parameter:
643c96a6f78SLois Curfman McInnes .  snes - SNES context
644c96a6f78SLois Curfman McInnes 
645c96a6f78SLois Curfman McInnes    Output Parameter:
646c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
647c96a6f78SLois Curfman McInnes 
648c96a6f78SLois Curfman McInnes    Notes:
649c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
650c96a6f78SLois Curfman McInnes 
65136851e7fSLois Curfman McInnes    Level: intermediate
65236851e7fSLois Curfman McInnes 
653c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
6542b668275SBarry Smith 
655*58ebbce7SBarry Smith .seealso:  SNESGetIterationNumber(), SNESGetFunctionNorm()S, NESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures()
656c96a6f78SLois Curfman McInnes @*/
657b850b91aSLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveIterations(SNES snes,PetscInt* lits)
658c96a6f78SLois Curfman McInnes {
6593a40ed3dSBarry Smith   PetscFunctionBegin;
6604482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6614482741eSBarry Smith   PetscValidIntPointer(lits,2);
662c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
664c96a6f78SLois Curfman McInnes }
665c96a6f78SLois Curfman McInnes 
6664a2ae208SSatish Balay #undef __FUNCT__
66794b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
66852baeb72SSatish Balay /*@
66994b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
6709b94acceSBarry Smith 
67194b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
672c7afd0dbSLois Curfman McInnes 
6739b94acceSBarry Smith    Input Parameter:
6749b94acceSBarry Smith .  snes - the SNES context
6759b94acceSBarry Smith 
6769b94acceSBarry Smith    Output Parameter:
67794b7f48cSBarry Smith .  ksp - the KSP context
6789b94acceSBarry Smith 
6799b94acceSBarry Smith    Notes:
68094b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
6819b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
6822999313aSBarry Smith    PC contexts as well.
6839b94acceSBarry Smith 
68436851e7fSLois Curfman McInnes    Level: beginner
68536851e7fSLois Curfman McInnes 
68694b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
6879b94acceSBarry Smith 
6882999313aSBarry Smith .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
6899b94acceSBarry Smith @*/
69063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp)
6919b94acceSBarry Smith {
6923a40ed3dSBarry Smith   PetscFunctionBegin;
6934482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6944482741eSBarry Smith   PetscValidPointer(ksp,2);
69594b7f48cSBarry Smith   *ksp = snes->ksp;
6963a40ed3dSBarry Smith   PetscFunctionReturn(0);
6979b94acceSBarry Smith }
69882bf6240SBarry Smith 
6994a2ae208SSatish Balay #undef __FUNCT__
7002999313aSBarry Smith #define __FUNCT__ "SNESSetKSP"
7012999313aSBarry Smith /*@
7022999313aSBarry Smith    SNESSetKSP - Sets a KSP context for the SNES object to use
7032999313aSBarry Smith 
7042999313aSBarry Smith    Not Collective, but the SNES and KSP objects must live on the same MPI_Comm
7052999313aSBarry Smith 
7062999313aSBarry Smith    Input Parameters:
7072999313aSBarry Smith +  snes - the SNES context
7082999313aSBarry Smith -  ksp - the KSP context
7092999313aSBarry Smith 
7102999313aSBarry Smith    Notes:
7112999313aSBarry Smith    The SNES object already has its KSP object, you can obtain with SNESGetKSP()
7122999313aSBarry Smith    so this routine is rarely needed.
7132999313aSBarry Smith 
7142999313aSBarry Smith    The KSP object that is already in the SNES object has its reference count
7152999313aSBarry Smith    decreased by one.
7162999313aSBarry Smith 
7172999313aSBarry Smith    Level: developer
7182999313aSBarry Smith 
7192999313aSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
7202999313aSBarry Smith 
7212999313aSBarry Smith .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP()
7222999313aSBarry Smith @*/
7232999313aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetKSP(SNES snes,KSP ksp)
7242999313aSBarry Smith {
7252999313aSBarry Smith   PetscErrorCode ierr;
7262999313aSBarry Smith 
7272999313aSBarry Smith   PetscFunctionBegin;
7282999313aSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7292999313aSBarry Smith   PetscValidHeaderSpecific(ksp,KSP_COOKIE,2);
7302999313aSBarry Smith   PetscCheckSameComm(snes,1,ksp,2);
7317dcf0eaaSdalcinl   ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr);
732906ed7ccSBarry Smith   if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);}
7332999313aSBarry Smith   snes->ksp = ksp;
7342999313aSBarry Smith   PetscFunctionReturn(0);
7352999313aSBarry Smith }
7362999313aSBarry Smith 
7377adad957SLisandro Dalcin #if 0
7382999313aSBarry Smith #undef __FUNCT__
7394a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
7406849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
741e24b481bSBarry Smith {
742e24b481bSBarry Smith   PetscFunctionBegin;
743e24b481bSBarry Smith   PetscFunctionReturn(0);
744e24b481bSBarry Smith }
7457adad957SLisandro Dalcin #endif
746e24b481bSBarry Smith 
7479b94acceSBarry Smith /* -----------------------------------------------------------*/
7484a2ae208SSatish Balay #undef __FUNCT__
7494a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
75052baeb72SSatish Balay /*@
7519b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
7529b94acceSBarry Smith 
753c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
754c7afd0dbSLois Curfman McInnes 
755c7afd0dbSLois Curfman McInnes    Input Parameters:
756906ed7ccSBarry Smith .  comm - MPI communicator
7579b94acceSBarry Smith 
7589b94acceSBarry Smith    Output Parameter:
7599b94acceSBarry Smith .  outsnes - the new SNES context
7609b94acceSBarry Smith 
761c7afd0dbSLois Curfman McInnes    Options Database Keys:
762c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
763c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
764c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
765c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
766c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
767c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
768c1f60f51SBarry Smith 
76936851e7fSLois Curfman McInnes    Level: beginner
77036851e7fSLois Curfman McInnes 
7719b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
7729b94acceSBarry Smith 
7734b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES
7749b94acceSBarry Smith @*/
77563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes)
7769b94acceSBarry Smith {
777dfbe8321SBarry Smith   PetscErrorCode      ierr;
7789b94acceSBarry Smith   SNES                snes;
779fa9f3622SBarry Smith   SNESKSPEW           *kctx;
78037fcc0dbSBarry Smith 
7813a40ed3dSBarry Smith   PetscFunctionBegin;
782ed1caa07SMatthew Knepley   PetscValidPointer(outsnes,2);
7838ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
7848ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
7858ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
7868ba1e511SMatthew Knepley #endif
7878ba1e511SMatthew Knepley 
788e7788613SBarry Smith   ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
7897adad957SLisandro Dalcin 
79085385478SLisandro Dalcin   snes->ops->converged    = SNESDefaultConverged;
7919b94acceSBarry Smith   snes->max_its           = 50;
7929750a799SBarry Smith   snes->max_funcs	  = 10000;
7939b94acceSBarry Smith   snes->norm		  = 0.0;
794b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
795b4874afaSBarry Smith   snes->ttol              = 0.0;
79670441072SBarry Smith   snes->abstol		  = 1.e-50;
7979b94acceSBarry Smith   snes->xtol		  = 1.e-8;
7984b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
7999b94acceSBarry Smith   snes->nfuncs            = 0;
80050ffb88aSMatthew Knepley   snes->numFailures       = 0;
80150ffb88aSMatthew Knepley   snes->maxFailures       = 1;
8027a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
803639f9d9dSBarry Smith   snes->numbermonitors    = 0;
8049b94acceSBarry Smith   snes->data              = 0;
8054dc4c822SBarry Smith   snes->setupcalled       = PETSC_FALSE;
806186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
8076f24a144SLois Curfman McInnes   snes->vwork             = 0;
8086f24a144SLois Curfman McInnes   snes->nwork             = 0;
809758f92a0SBarry Smith   snes->conv_hist_len     = 0;
810758f92a0SBarry Smith   snes->conv_hist_max     = 0;
811758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
812758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
813758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
814184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
8159b94acceSBarry Smith 
8163d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
8173d4c4710SBarry Smith   snes->maxLinearSolveFailures = 1;
8183d4c4710SBarry Smith 
8199b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
82038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr);
8219b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
8229b94acceSBarry Smith   kctx->version     = 2;
8239b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
8249b94acceSBarry Smith                              this was too large for some test cases */
8259b94acceSBarry Smith   kctx->rtol_last   = 0;
8269b94acceSBarry Smith   kctx->rtol_max    = .9;
8279b94acceSBarry Smith   kctx->gamma       = 1.0;
82871f87433Sdalcinl   kctx->alpha       = .5*(1.0 + sqrt(5.0));
82971f87433Sdalcinl   kctx->alpha2      = kctx->alpha;
8309b94acceSBarry Smith   kctx->threshold   = .1;
8319b94acceSBarry Smith   kctx->lresid_last = 0;
8329b94acceSBarry Smith   kctx->norm_last   = 0;
8339b94acceSBarry Smith 
83494b7f48cSBarry Smith   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
83552e6d16bSBarry Smith   ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
8365334005bSBarry Smith 
8379b94acceSBarry Smith   *outsnes = snes;
83800036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
8393a40ed3dSBarry Smith   PetscFunctionReturn(0);
8409b94acceSBarry Smith }
8419b94acceSBarry Smith 
8424a2ae208SSatish Balay #undef __FUNCT__
8434a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
8449b94acceSBarry Smith /*@C
8459b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
8469b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
8479b94acceSBarry Smith    equations.
8489b94acceSBarry Smith 
849fee21e36SBarry Smith    Collective on SNES
850fee21e36SBarry Smith 
851c7afd0dbSLois Curfman McInnes    Input Parameters:
852c7afd0dbSLois Curfman McInnes +  snes - the SNES context
853c7afd0dbSLois Curfman McInnes .  r - vector to store function value
854de044059SHong Zhang .  func - function evaluation routine
855c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
856c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
8579b94acceSBarry Smith 
858c7afd0dbSLois Curfman McInnes    Calling sequence of func:
8598d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
860c7afd0dbSLois Curfman McInnes 
861313e4042SLois Curfman McInnes .  f - function vector
862c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
8639b94acceSBarry Smith 
8649b94acceSBarry Smith    Notes:
8659b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
8669b94acceSBarry Smith $      f'(x) x = -f(x),
867c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
8689b94acceSBarry Smith 
86936851e7fSLois Curfman McInnes    Level: beginner
87036851e7fSLois Curfman McInnes 
8719b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8729b94acceSBarry Smith 
873a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
8749b94acceSBarry Smith @*/
87563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
8769b94acceSBarry Smith {
87785385478SLisandro Dalcin   PetscErrorCode ierr;
8783a40ed3dSBarry Smith   PetscFunctionBegin;
8794482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8804482741eSBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
881c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,r,2);
88285385478SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr);
88385385478SLisandro Dalcin   if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); }
884e7788613SBarry Smith   snes->ops->computefunction = func;
88585385478SLisandro Dalcin   snes->vec_func             = r;
8869b94acceSBarry Smith   snes->funP                 = ctx;
8873a40ed3dSBarry Smith   PetscFunctionReturn(0);
8889b94acceSBarry Smith }
8899b94acceSBarry Smith 
8903ab0aad5SBarry Smith /* --------------------------------------------------------------- */
8913ab0aad5SBarry Smith #undef __FUNCT__
8921096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs"
8931096aae1SMatthew Knepley /*@C
8941096aae1SMatthew Knepley    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
8951096aae1SMatthew Knepley    it assumes a zero right hand side.
8961096aae1SMatthew Knepley 
8971096aae1SMatthew Knepley    Collective on SNES
8981096aae1SMatthew Knepley 
8991096aae1SMatthew Knepley    Input Parameter:
9001096aae1SMatthew Knepley .  snes - the SNES context
9011096aae1SMatthew Knepley 
9021096aae1SMatthew Knepley    Output Parameter:
903bc08b0f1SBarry Smith .  rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null
9041096aae1SMatthew Knepley 
9051096aae1SMatthew Knepley    Level: intermediate
9061096aae1SMatthew Knepley 
9071096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side
9081096aae1SMatthew Knepley 
90985385478SLisandro Dalcin .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
9101096aae1SMatthew Knepley @*/
9111096aae1SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs)
9121096aae1SMatthew Knepley {
9131096aae1SMatthew Knepley   PetscFunctionBegin;
9141096aae1SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9151096aae1SMatthew Knepley   PetscValidPointer(rhs,2);
91685385478SLisandro Dalcin   *rhs = snes->vec_rhs;
9171096aae1SMatthew Knepley   PetscFunctionReturn(0);
9181096aae1SMatthew Knepley }
9191096aae1SMatthew Knepley 
9201096aae1SMatthew Knepley #undef __FUNCT__
9214a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
9229b94acceSBarry Smith /*@
92336851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
9249b94acceSBarry Smith                          SNESSetFunction().
9259b94acceSBarry Smith 
926c7afd0dbSLois Curfman McInnes    Collective on SNES
927c7afd0dbSLois Curfman McInnes 
9289b94acceSBarry Smith    Input Parameters:
929c7afd0dbSLois Curfman McInnes +  snes - the SNES context
930c7afd0dbSLois Curfman McInnes -  x - input vector
9319b94acceSBarry Smith 
9329b94acceSBarry Smith    Output Parameter:
9333638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
9349b94acceSBarry Smith 
9351bffabb2SLois Curfman McInnes    Notes:
93636851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
93736851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
93836851e7fSLois Curfman McInnes    themselves.
93936851e7fSLois Curfman McInnes 
94036851e7fSLois Curfman McInnes    Level: developer
94136851e7fSLois Curfman McInnes 
9429b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
9439b94acceSBarry Smith 
944a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
9459b94acceSBarry Smith @*/
94663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y)
9479b94acceSBarry Smith {
948dfbe8321SBarry Smith   PetscErrorCode ierr;
9499b94acceSBarry Smith 
9503a40ed3dSBarry Smith   PetscFunctionBegin;
9514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9524482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
9534482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
954c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
955c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
956184914b5SBarry Smith 
957d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
958e7788613SBarry Smith   if (snes->ops->computefunction) {
959d64ed03dSBarry Smith     PetscStackPush("SNES user function");
960e9a2bbcdSBarry Smith     CHKMEMQ;
961e7788613SBarry Smith     ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP);
962e9a2bbcdSBarry Smith     CHKMEMQ;
963d64ed03dSBarry Smith     PetscStackPop;
964d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
96519717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
96619717074SBarry Smith     }
967d5e45103SBarry Smith     CHKERRQ(ierr);
96885385478SLisandro Dalcin   } else if (snes->vec_rhs) {
9691096aae1SMatthew Knepley     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
9701096aae1SMatthew Knepley   } else {
9711096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
9721096aae1SMatthew Knepley   }
97385385478SLisandro Dalcin   if (snes->vec_rhs) {
97485385478SLisandro Dalcin     ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr);
9753ab0aad5SBarry Smith   }
976ae3c334cSLois Curfman McInnes   snes->nfuncs++;
977d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
9783a40ed3dSBarry Smith   PetscFunctionReturn(0);
9799b94acceSBarry Smith }
9809b94acceSBarry Smith 
9814a2ae208SSatish Balay #undef __FUNCT__
9824a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
98362fef451SLois Curfman McInnes /*@
98462fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
98562fef451SLois Curfman McInnes    set with SNESSetJacobian().
98662fef451SLois Curfman McInnes 
987c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
988c7afd0dbSLois Curfman McInnes 
98962fef451SLois Curfman McInnes    Input Parameters:
990c7afd0dbSLois Curfman McInnes +  snes - the SNES context
991c7afd0dbSLois Curfman McInnes -  x - input vector
99262fef451SLois Curfman McInnes 
99362fef451SLois Curfman McInnes    Output Parameters:
994c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
99562fef451SLois Curfman McInnes .  B - optional preconditioning matrix
9962b668275SBarry Smith -  flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER)
997fee21e36SBarry Smith 
99862fef451SLois Curfman McInnes    Notes:
99962fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
100062fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
100162fef451SLois Curfman McInnes 
100294b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
1003dc5a77f8SLois Curfman McInnes    flag parameter.
100462fef451SLois Curfman McInnes 
100536851e7fSLois Curfman McInnes    Level: developer
100636851e7fSLois Curfman McInnes 
100762fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
100862fef451SLois Curfman McInnes 
10092b668275SBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators(), MatStructure
101062fef451SLois Curfman McInnes @*/
101163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
10129b94acceSBarry Smith {
1013dfbe8321SBarry Smith   PetscErrorCode ierr;
10143a40ed3dSBarry Smith 
10153a40ed3dSBarry Smith   PetscFunctionBegin;
10164482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10174482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
10184482741eSBarry Smith   PetscValidPointer(flg,5);
1019c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
1020e7788613SBarry Smith   if (!snes->ops->computejacobian) PetscFunctionReturn(0);
1021d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1022c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
1023d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
1024dc67937bSBarry Smith   CHKMEMQ;
1025e7788613SBarry Smith   ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1026dc67937bSBarry Smith   CHKMEMQ;
1027d64ed03dSBarry Smith   PetscStackPop;
1028d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
10296d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
10306ce558aeSBarry Smith   /* PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
10316ce558aeSBarry Smith     PetscValidHeaderSpecific(*B,MAT_COOKIE,4);   */
10323a40ed3dSBarry Smith   PetscFunctionReturn(0);
10339b94acceSBarry Smith }
10349b94acceSBarry Smith 
10354a2ae208SSatish Balay #undef __FUNCT__
10364a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
10379b94acceSBarry Smith /*@C
10389b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1039044dda88SLois Curfman McInnes    location to store the matrix.
10409b94acceSBarry Smith 
1041c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1042c7afd0dbSLois Curfman McInnes 
10439b94acceSBarry Smith    Input Parameters:
1044c7afd0dbSLois Curfman McInnes +  snes - the SNES context
10459b94acceSBarry Smith .  A - Jacobian matrix
10469b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
10479b94acceSBarry Smith .  func - Jacobian evaluation routine
1048c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
10492cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
10509b94acceSBarry Smith 
10519b94acceSBarry Smith    Calling sequence of func:
10528d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10539b94acceSBarry Smith 
1054c7afd0dbSLois Curfman McInnes +  x - input vector
10559b94acceSBarry Smith .  A - Jacobian matrix
10569b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1057ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
10582b668275SBarry Smith    structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
1059c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
10609b94acceSBarry Smith 
10619b94acceSBarry Smith    Notes:
106294b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
10632cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1064ac21db08SLois Curfman McInnes 
1065ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
10669b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
10679b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10689b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10699b94acceSBarry Smith    throughout the global iterations.
10709b94acceSBarry Smith 
107136851e7fSLois Curfman McInnes    Level: beginner
107236851e7fSLois Curfman McInnes 
10739b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
10749b94acceSBarry Smith 
10753ec795f1SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure
10769b94acceSBarry Smith @*/
107763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
10789b94acceSBarry Smith {
1079dfbe8321SBarry Smith   PetscErrorCode ierr;
10803a7fca6bSBarry Smith 
10813a40ed3dSBarry Smith   PetscFunctionBegin;
10824482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10834482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
10844482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
1085c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
1086c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
1087e7788613SBarry Smith   if (func) snes->ops->computejacobian = func;
10883a7fca6bSBarry Smith   if (ctx)  snes->jacP                 = ctx;
10893a7fca6bSBarry Smith   if (A) {
10907dcf0eaaSdalcinl     ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
10913a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
10929b94acceSBarry Smith     snes->jacobian = A;
10933a7fca6bSBarry Smith   }
10943a7fca6bSBarry Smith   if (B) {
10957dcf0eaaSdalcinl     ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
10963a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
10979b94acceSBarry Smith     snes->jacobian_pre = B;
10983a7fca6bSBarry Smith   }
10993a40ed3dSBarry Smith   PetscFunctionReturn(0);
11009b94acceSBarry Smith }
110162fef451SLois Curfman McInnes 
11024a2ae208SSatish Balay #undef __FUNCT__
11034a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
1104c2aafc4cSSatish Balay /*@C
1105b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1106b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1107b4fd4287SBarry Smith 
1108c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1109c7afd0dbSLois Curfman McInnes 
1110b4fd4287SBarry Smith    Input Parameter:
1111b4fd4287SBarry Smith .  snes - the nonlinear solver context
1112b4fd4287SBarry Smith 
1113b4fd4287SBarry Smith    Output Parameters:
1114c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1115b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
111670e92668SMatthew Knepley .  func - location to put Jacobian function (or PETSC_NULL)
111770e92668SMatthew Knepley -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
1118fee21e36SBarry Smith 
111936851e7fSLois Curfman McInnes    Level: advanced
112036851e7fSLois Curfman McInnes 
1121b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1122b4fd4287SBarry Smith @*/
112370e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
1124b4fd4287SBarry Smith {
11253a40ed3dSBarry Smith   PetscFunctionBegin;
11264482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1127b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
1128b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
1129e7788613SBarry Smith   if (func) *func = snes->ops->computejacobian;
113070e92668SMatthew Knepley   if (ctx)  *ctx  = snes->jacP;
11313a40ed3dSBarry Smith   PetscFunctionReturn(0);
1132b4fd4287SBarry Smith }
1133b4fd4287SBarry Smith 
11349b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
113563dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
11369b94acceSBarry Smith 
11374a2ae208SSatish Balay #undef __FUNCT__
11384a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
11399b94acceSBarry Smith /*@
11409b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1141272ac6f2SLois Curfman McInnes    of a nonlinear solver.
11429b94acceSBarry Smith 
1143fee21e36SBarry Smith    Collective on SNES
1144fee21e36SBarry Smith 
1145c7afd0dbSLois Curfman McInnes    Input Parameters:
114670e92668SMatthew Knepley .  snes - the SNES context
1147c7afd0dbSLois Curfman McInnes 
1148272ac6f2SLois Curfman McInnes    Notes:
1149272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1150272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1151272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1152272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1153272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1154272ac6f2SLois Curfman McInnes 
115536851e7fSLois Curfman McInnes    Level: advanced
115636851e7fSLois Curfman McInnes 
11579b94acceSBarry Smith .keywords: SNES, nonlinear, setup
11589b94acceSBarry Smith 
11599b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
11609b94acceSBarry Smith @*/
116170e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
11629b94acceSBarry Smith {
1163dfbe8321SBarry Smith   PetscErrorCode ierr;
116471f87433Sdalcinl   PetscTruth     flg;
11653a40ed3dSBarry Smith 
11663a40ed3dSBarry Smith   PetscFunctionBegin;
11674482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
11684dc4c822SBarry Smith   if (snes->setupcalled) PetscFunctionReturn(0);
11699b94acceSBarry Smith 
11707adad957SLisandro Dalcin   if (!((PetscObject)snes)->type_name) {
117185385478SLisandro Dalcin     ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr);
117285385478SLisandro Dalcin   }
117385385478SLisandro Dalcin 
11747adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1175c1f60f51SBarry Smith   /*
1176c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1177dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1178c1f60f51SBarry Smith   */
1179c1f60f51SBarry Smith   if (flg) {
1180c1f60f51SBarry Smith     Mat J;
1181fef1beadSBarry Smith     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
11823ec795f1SBarry Smith     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
1183ae15b995SBarry Smith     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
11843a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
11853a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1186c1f60f51SBarry Smith   }
118745fc7adcSBarry Smith 
118803c60df9SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT)
11897adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
119045fc7adcSBarry Smith   if (flg) {
119145fc7adcSBarry Smith     Mat J;
119245fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
119375396ef9SLisandro Dalcin     ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr);
119445fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
119545fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
119645fc7adcSBarry Smith   }
119732a4b47aSMatthew Knepley #endif
119845fc7adcSBarry Smith 
11997adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1200c1f60f51SBarry Smith   /*
1201dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1202c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1203c1f60f51SBarry Smith    */
120431615d04SLois Curfman McInnes   if (flg) {
1205272ac6f2SLois Curfman McInnes     Mat  J;
1206b5d62d44SBarry Smith     KSP ksp;
120794b7f48cSBarry Smith     PC   pc;
120875396ef9SLisandro Dalcin     /* create and set matrix-free operator */
1209fef1beadSBarry Smith     ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr);
12103ec795f1SBarry Smith     ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr);
121175396ef9SLisandro Dalcin     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
12123ec795f1SBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr);
12133a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1214f3ef73ceSBarry Smith     /* force no preconditioner */
121594b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1216b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1217a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1218a9815358SBarry Smith     if (!flg) {
121975396ef9SLisandro Dalcin       ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr);
1220f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1221272ac6f2SLois Curfman McInnes     }
1222a9815358SBarry Smith   }
1223f3ef73ceSBarry Smith 
122485385478SLisandro Dalcin   if (!snes->vec_func && !snes->vec_rhs) {
12251096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
12261096aae1SMatthew Knepley   }
122785385478SLisandro Dalcin   if (!snes->ops->computefunction && !snes->vec_rhs) {
12281096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
12291096aae1SMatthew Knepley   }
123075396ef9SLisandro Dalcin   if (!snes->jacobian) {
123175396ef9SLisandro Dalcin     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
123275396ef9SLisandro Dalcin   }
1233a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
123429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1235a8c6a408SBarry Smith   }
1236a703fe33SLois Curfman McInnes 
1237410397dcSLisandro Dalcin   if (snes->ops->setup) {
1238410397dcSLisandro Dalcin     ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr);
1239410397dcSLisandro Dalcin   }
12407aaed0d8SBarry Smith   snes->setupcalled = PETSC_TRUE;
12413a40ed3dSBarry Smith   PetscFunctionReturn(0);
12429b94acceSBarry Smith }
12439b94acceSBarry Smith 
12444a2ae208SSatish Balay #undef __FUNCT__
12454a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
124652baeb72SSatish Balay /*@
12479b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
12489b94acceSBarry Smith    with SNESCreate().
12499b94acceSBarry Smith 
1250c7afd0dbSLois Curfman McInnes    Collective on SNES
1251c7afd0dbSLois Curfman McInnes 
12529b94acceSBarry Smith    Input Parameter:
12539b94acceSBarry Smith .  snes - the SNES context
12549b94acceSBarry Smith 
125536851e7fSLois Curfman McInnes    Level: beginner
125636851e7fSLois Curfman McInnes 
12579b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
12589b94acceSBarry Smith 
125963a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
12609b94acceSBarry Smith @*/
126163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
12629b94acceSBarry Smith {
12636849ba73SBarry Smith   PetscErrorCode ierr;
12643a40ed3dSBarry Smith 
12653a40ed3dSBarry Smith   PetscFunctionBegin;
12664482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12677adad957SLisandro Dalcin   if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0);
1268d4bb536fSBarry Smith 
1269be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
12700f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1271e7788613SBarry Smith   if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);}
127285385478SLisandro Dalcin 
127385385478SLisandro Dalcin   if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);}
127485385478SLisandro Dalcin   if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);}
127585385478SLisandro Dalcin   if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);}
12763a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
12773a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
127894b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
127985385478SLisandro Dalcin   ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);
1280522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1281a6570f20SBarry Smith   ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
1282a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
12833a40ed3dSBarry Smith   PetscFunctionReturn(0);
12849b94acceSBarry Smith }
12859b94acceSBarry Smith 
12869b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
12879b94acceSBarry Smith 
12884a2ae208SSatish Balay #undef __FUNCT__
12894a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
12909b94acceSBarry Smith /*@
1291d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
12929b94acceSBarry Smith 
1293c7afd0dbSLois Curfman McInnes    Collective on SNES
1294c7afd0dbSLois Curfman McInnes 
12959b94acceSBarry Smith    Input Parameters:
1296c7afd0dbSLois Curfman McInnes +  snes - the SNES context
129770441072SBarry Smith .  abstol - absolute convergence tolerance
129833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
129933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
130033174efeSLois Curfman McInnes            of the change in the solution between steps
130133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1302c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1303fee21e36SBarry Smith 
130433174efeSLois Curfman McInnes    Options Database Keys:
130570441072SBarry Smith +    -snes_atol <abstol> - Sets abstol
1306c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1307c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1308c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1309c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
13109b94acceSBarry Smith 
1311d7a720efSLois Curfman McInnes    Notes:
13129b94acceSBarry Smith    The default maximum number of iterations is 50.
13139b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
13149b94acceSBarry Smith 
131536851e7fSLois Curfman McInnes    Level: intermediate
131636851e7fSLois Curfman McInnes 
131733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
13189b94acceSBarry Smith 
13192492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance()
13209b94acceSBarry Smith @*/
132163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
13229b94acceSBarry Smith {
13233a40ed3dSBarry Smith   PetscFunctionBegin;
13244482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
132570441072SBarry Smith   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1326d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1327d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1328d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1329d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
13303a40ed3dSBarry Smith   PetscFunctionReturn(0);
13319b94acceSBarry Smith }
13329b94acceSBarry Smith 
13334a2ae208SSatish Balay #undef __FUNCT__
13344a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
13359b94acceSBarry Smith /*@
133633174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
133733174efeSLois Curfman McInnes 
1338c7afd0dbSLois Curfman McInnes    Not Collective
1339c7afd0dbSLois Curfman McInnes 
134033174efeSLois Curfman McInnes    Input Parameters:
1341c7afd0dbSLois Curfman McInnes +  snes - the SNES context
134285385478SLisandro Dalcin .  atol - absolute convergence tolerance
134333174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
134433174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
134533174efeSLois Curfman McInnes            of the change in the solution between steps
134633174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1347c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1348fee21e36SBarry Smith 
134933174efeSLois Curfman McInnes    Notes:
135033174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
135133174efeSLois Curfman McInnes 
135236851e7fSLois Curfman McInnes    Level: intermediate
135336851e7fSLois Curfman McInnes 
135433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
135533174efeSLois Curfman McInnes 
135633174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
135733174efeSLois Curfman McInnes @*/
135885385478SLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
135933174efeSLois Curfman McInnes {
13603a40ed3dSBarry Smith   PetscFunctionBegin;
13614482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
136285385478SLisandro Dalcin   if (atol)  *atol  = snes->abstol;
136333174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
136433174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
136533174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
136633174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
13673a40ed3dSBarry Smith   PetscFunctionReturn(0);
136833174efeSLois Curfman McInnes }
136933174efeSLois Curfman McInnes 
13704a2ae208SSatish Balay #undef __FUNCT__
13714a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
137233174efeSLois Curfman McInnes /*@
13739b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
13749b94acceSBarry Smith 
1375fee21e36SBarry Smith    Collective on SNES
1376fee21e36SBarry Smith 
1377c7afd0dbSLois Curfman McInnes    Input Parameters:
1378c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1379c7afd0dbSLois Curfman McInnes -  tol - tolerance
1380c7afd0dbSLois Curfman McInnes 
13819b94acceSBarry Smith    Options Database Key:
1382c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
13839b94acceSBarry Smith 
138436851e7fSLois Curfman McInnes    Level: intermediate
138536851e7fSLois Curfman McInnes 
13869b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
13879b94acceSBarry Smith 
13882492ecdbSBarry Smith .seealso: SNESSetTolerances()
13899b94acceSBarry Smith @*/
139063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
13919b94acceSBarry Smith {
13923a40ed3dSBarry Smith   PetscFunctionBegin;
13934482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13949b94acceSBarry Smith   snes->deltatol = tol;
13953a40ed3dSBarry Smith   PetscFunctionReturn(0);
13969b94acceSBarry Smith }
13979b94acceSBarry Smith 
1398df9fa365SBarry Smith /*
1399df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1400df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1401df9fa365SBarry Smith    macros instead of functions
1402df9fa365SBarry Smith */
14034a2ae208SSatish Balay #undef __FUNCT__
1404a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLG"
1405a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1406ce1608b8SBarry Smith {
1407dfbe8321SBarry Smith   PetscErrorCode ierr;
1408ce1608b8SBarry Smith 
1409ce1608b8SBarry Smith   PetscFunctionBegin;
14104482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1411a6570f20SBarry Smith   ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1412ce1608b8SBarry Smith   PetscFunctionReturn(0);
1413ce1608b8SBarry Smith }
1414ce1608b8SBarry Smith 
14154a2ae208SSatish Balay #undef __FUNCT__
1416a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLGCreate"
1417a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1418df9fa365SBarry Smith {
1419dfbe8321SBarry Smith   PetscErrorCode ierr;
1420df9fa365SBarry Smith 
1421df9fa365SBarry Smith   PetscFunctionBegin;
1422a6570f20SBarry Smith   ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1423df9fa365SBarry Smith   PetscFunctionReturn(0);
1424df9fa365SBarry Smith }
1425df9fa365SBarry Smith 
14264a2ae208SSatish Balay #undef __FUNCT__
1427a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorLGDestroy"
1428a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw)
1429df9fa365SBarry Smith {
1430dfbe8321SBarry Smith   PetscErrorCode ierr;
1431df9fa365SBarry Smith 
1432df9fa365SBarry Smith   PetscFunctionBegin;
1433a6570f20SBarry Smith   ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr);
1434df9fa365SBarry Smith   PetscFunctionReturn(0);
1435df9fa365SBarry Smith }
1436df9fa365SBarry Smith 
14379b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
14389b94acceSBarry Smith 
14394a2ae208SSatish Balay #undef __FUNCT__
1440a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorSet"
14419b94acceSBarry Smith /*@C
1442a6570f20SBarry Smith    SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every
14439b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
14449b94acceSBarry Smith    progress.
14459b94acceSBarry Smith 
1446fee21e36SBarry Smith    Collective on SNES
1447fee21e36SBarry Smith 
1448c7afd0dbSLois Curfman McInnes    Input Parameters:
1449c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1450c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1451b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1452e8105e01SRichard Katz           monitor routine (use PETSC_NULL if no context is desired)
1453b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1454b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
14559b94acceSBarry Smith 
1456c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1457a7cc72afSBarry Smith $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1458c7afd0dbSLois Curfman McInnes 
1459c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1460c7afd0dbSLois Curfman McInnes .    its - iteration number
1461c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
146240a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
14639b94acceSBarry Smith 
14649665c990SLois Curfman McInnes    Options Database Keys:
1465a6570f20SBarry Smith +    -snes_monitor        - sets SNESMonitorDefault()
1466a6570f20SBarry Smith .    -snes_monitor_draw    - sets line graph monitor,
1467a6570f20SBarry Smith                             uses SNESMonitorLGCreate()
1468a6570f20SBarry Smith _    -snes_monitor_cancel - cancels all monitors that have
1469c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1470a6570f20SBarry Smith                             calls to SNESMonitorSet(), but
1471c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1472c7afd0dbSLois Curfman McInnes                             the options database.
14739665c990SLois Curfman McInnes 
1474639f9d9dSBarry Smith    Notes:
14756bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
1476a6570f20SBarry Smith    SNESMonitorSet() multiple times; all will be called in the
14776bc08f3fSLois Curfman McInnes    order in which they were set.
1478639f9d9dSBarry Smith 
147936851e7fSLois Curfman McInnes    Level: intermediate
148036851e7fSLois Curfman McInnes 
14819b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
14829b94acceSBarry Smith 
1483a6570f20SBarry Smith .seealso: SNESMonitorDefault(), SNESMonitorCancel()
14849b94acceSBarry Smith @*/
1485b90d0a6eSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
14869b94acceSBarry Smith {
1487b90d0a6eSBarry Smith   PetscInt i;
1488b90d0a6eSBarry Smith 
14893a40ed3dSBarry Smith   PetscFunctionBegin;
14904482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1491639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
149229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1493639f9d9dSBarry Smith   }
1494b90d0a6eSBarry Smith   for (i=0; i<snes->numbermonitors;i++) {
1495b90d0a6eSBarry Smith     if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0);
1496b90d0a6eSBarry Smith 
1497b90d0a6eSBarry Smith     /* check if both default monitors that share common ASCII viewer */
1498b90d0a6eSBarry Smith     if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) {
1499b90d0a6eSBarry Smith       if (mctx && snes->monitorcontext[i]) {
1500b90d0a6eSBarry Smith         PetscErrorCode          ierr;
1501b90d0a6eSBarry Smith         PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx;
1502b90d0a6eSBarry Smith         PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i];
1503b90d0a6eSBarry Smith         if (viewer1->viewer == viewer2->viewer) {
1504b90d0a6eSBarry Smith           ierr = (*monitordestroy)(mctx);CHKERRQ(ierr);
1505b90d0a6eSBarry Smith           PetscFunctionReturn(0);
1506b90d0a6eSBarry Smith         }
1507b90d0a6eSBarry Smith       }
1508b90d0a6eSBarry Smith     }
1509b90d0a6eSBarry Smith   }
1510b90d0a6eSBarry Smith   snes->monitor[snes->numbermonitors]           = monitor;
1511b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1512639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
15133a40ed3dSBarry Smith   PetscFunctionReturn(0);
15149b94acceSBarry Smith }
15159b94acceSBarry Smith 
15164a2ae208SSatish Balay #undef __FUNCT__
1517a6570f20SBarry Smith #define __FUNCT__ "SNESMonitorCancel"
15185cd90555SBarry Smith /*@C
1519a6570f20SBarry Smith    SNESMonitorCancel - Clears all the monitor functions for a SNES object.
15205cd90555SBarry Smith 
1521c7afd0dbSLois Curfman McInnes    Collective on SNES
1522c7afd0dbSLois Curfman McInnes 
15235cd90555SBarry Smith    Input Parameters:
15245cd90555SBarry Smith .  snes - the SNES context
15255cd90555SBarry Smith 
15261a480d89SAdministrator    Options Database Key:
1527a6570f20SBarry Smith .  -snes_monitor_cancel - cancels all monitors that have been hardwired
1528a6570f20SBarry Smith     into a code by calls to SNESMonitorSet(), but does not cancel those
1529c7afd0dbSLois Curfman McInnes     set via the options database
15305cd90555SBarry Smith 
15315cd90555SBarry Smith    Notes:
15325cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
15335cd90555SBarry Smith 
153436851e7fSLois Curfman McInnes    Level: intermediate
153536851e7fSLois Curfman McInnes 
15365cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
15375cd90555SBarry Smith 
1538a6570f20SBarry Smith .seealso: SNESMonitorDefault(), SNESMonitorSet()
15395cd90555SBarry Smith @*/
1540a6570f20SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes)
15415cd90555SBarry Smith {
1542d952e501SBarry Smith   PetscErrorCode ierr;
1543d952e501SBarry Smith   PetscInt       i;
1544d952e501SBarry Smith 
15455cd90555SBarry Smith   PetscFunctionBegin;
15464482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1547d952e501SBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1548d952e501SBarry Smith     if (snes->monitordestroy[i]) {
1549d952e501SBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1550d952e501SBarry Smith     }
1551d952e501SBarry Smith   }
15525cd90555SBarry Smith   snes->numbermonitors = 0;
15535cd90555SBarry Smith   PetscFunctionReturn(0);
15545cd90555SBarry Smith }
15555cd90555SBarry Smith 
15564a2ae208SSatish Balay #undef __FUNCT__
15574a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
15589b94acceSBarry Smith /*@C
15599b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
15609b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
15619b94acceSBarry Smith 
1562fee21e36SBarry Smith    Collective on SNES
1563fee21e36SBarry Smith 
1564c7afd0dbSLois Curfman McInnes    Input Parameters:
1565c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1566c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1567c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1568c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
15699b94acceSBarry Smith 
1570c7afd0dbSLois Curfman McInnes    Calling sequence of func:
157106ee9f85SBarry Smith $     PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1572c7afd0dbSLois Curfman McInnes 
1573c7afd0dbSLois Curfman McInnes +    snes - the SNES context
157406ee9f85SBarry Smith .    it - current iteration (0 is the first and is before any Newton step)
1575c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1576184914b5SBarry Smith .    reason - reason for convergence/divergence
1577c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
15784b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
15794b27c08aSLois Curfman McInnes -    f - 2-norm of function
15809b94acceSBarry Smith 
158136851e7fSLois Curfman McInnes    Level: advanced
158236851e7fSLois Curfman McInnes 
15839b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
15849b94acceSBarry Smith 
158585385478SLisandro Dalcin .seealso: SNESDefaultConverged(), SNESSkipConverged()
15869b94acceSBarry Smith @*/
158706ee9f85SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
15889b94acceSBarry Smith {
15893a40ed3dSBarry Smith   PetscFunctionBegin;
15904482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
159185385478SLisandro Dalcin   if (!func) func = SNESSkipConverged;
159285385478SLisandro Dalcin   snes->ops->converged = func;
159385385478SLisandro Dalcin   snes->cnvP           = cctx;
15943a40ed3dSBarry Smith   PetscFunctionReturn(0);
15959b94acceSBarry Smith }
15969b94acceSBarry Smith 
15974a2ae208SSatish Balay #undef __FUNCT__
15984a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
159952baeb72SSatish Balay /*@
1600184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1601184914b5SBarry Smith 
1602184914b5SBarry Smith    Not Collective
1603184914b5SBarry Smith 
1604184914b5SBarry Smith    Input Parameter:
1605184914b5SBarry Smith .  snes - the SNES context
1606184914b5SBarry Smith 
1607184914b5SBarry Smith    Output Parameter:
1608e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1609184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1610184914b5SBarry Smith 
1611184914b5SBarry Smith    Level: intermediate
1612184914b5SBarry Smith 
1613184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1614184914b5SBarry Smith 
1615184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1616184914b5SBarry Smith 
161785385478SLisandro Dalcin .seealso: SNESSetConvergenceTest(), SNESConvergedReason
1618184914b5SBarry Smith @*/
161963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1620184914b5SBarry Smith {
1621184914b5SBarry Smith   PetscFunctionBegin;
16224482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16234482741eSBarry Smith   PetscValidPointer(reason,2);
1624184914b5SBarry Smith   *reason = snes->reason;
1625184914b5SBarry Smith   PetscFunctionReturn(0);
1626184914b5SBarry Smith }
1627184914b5SBarry Smith 
16284a2ae208SSatish Balay #undef __FUNCT__
16294a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1630c9005455SLois Curfman McInnes /*@
1631c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1632c9005455SLois Curfman McInnes 
1633fee21e36SBarry Smith    Collective on SNES
1634fee21e36SBarry Smith 
1635c7afd0dbSLois Curfman McInnes    Input Parameters:
1636c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1637c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1638cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1639758f92a0SBarry Smith .  na  - size of a and its
164064731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1641758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1642c7afd0dbSLois Curfman McInnes 
1643c9005455SLois Curfman McInnes    Notes:
16444b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1645c9005455SLois Curfman McInnes    at each step.
1646c9005455SLois Curfman McInnes 
1647c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1648c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1649c9005455SLois Curfman McInnes    during the section of code that is being timed.
1650c9005455SLois Curfman McInnes 
165136851e7fSLois Curfman McInnes    Level: intermediate
165236851e7fSLois Curfman McInnes 
1653c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1654758f92a0SBarry Smith 
165508405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1656758f92a0SBarry Smith 
1657c9005455SLois Curfman McInnes @*/
1658a562a398SLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset)
1659c9005455SLois Curfman McInnes {
16603a40ed3dSBarry Smith   PetscFunctionBegin;
16614482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16624482741eSBarry Smith   if (na)  PetscValidScalarPointer(a,2);
1663a562a398SLisandro Dalcin   if (its) PetscValidIntPointer(its,3);
1664c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1665758f92a0SBarry Smith   snes->conv_hist_its   = its;
1666758f92a0SBarry Smith   snes->conv_hist_max   = na;
1667758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1668758f92a0SBarry Smith   PetscFunctionReturn(0);
1669758f92a0SBarry Smith }
1670758f92a0SBarry Smith 
16714a2ae208SSatish Balay #undef __FUNCT__
16724a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
16730c4c9dddSBarry Smith /*@C
1674758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1675758f92a0SBarry Smith 
1676758f92a0SBarry Smith    Collective on SNES
1677758f92a0SBarry Smith 
1678758f92a0SBarry Smith    Input Parameter:
1679758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1680758f92a0SBarry Smith 
1681758f92a0SBarry Smith    Output Parameters:
1682758f92a0SBarry Smith .  a   - array to hold history
1683758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1684758f92a0SBarry Smith          negative if not converged) for each solve.
1685758f92a0SBarry Smith -  na  - size of a and its
1686758f92a0SBarry Smith 
1687758f92a0SBarry Smith    Notes:
1688758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1689758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1690758f92a0SBarry Smith 
1691758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1692758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1693758f92a0SBarry Smith    during the section of code that is being timed.
1694758f92a0SBarry Smith 
1695758f92a0SBarry Smith    Level: intermediate
1696758f92a0SBarry Smith 
1697758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1698758f92a0SBarry Smith 
1699758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1700758f92a0SBarry Smith 
1701758f92a0SBarry Smith @*/
170263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1703758f92a0SBarry Smith {
1704758f92a0SBarry Smith   PetscFunctionBegin;
17054482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1706758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1707758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1708758f92a0SBarry Smith   if (na)  *na  = snes->conv_hist_len;
17093a40ed3dSBarry Smith   PetscFunctionReturn(0);
1710c9005455SLois Curfman McInnes }
1711c9005455SLois Curfman McInnes 
1712e74ef692SMatthew Knepley #undef __FUNCT__
1713e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1714ac226902SBarry Smith /*@C
171576b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
17167e4bb74cSBarry Smith   at the beginning o every iteration of the nonlinear solve. Specifically
17177e4bb74cSBarry Smith   it is called just before the Jacobian is "evaluated".
171876b2cf59SMatthew Knepley 
171976b2cf59SMatthew Knepley   Collective on SNES
172076b2cf59SMatthew Knepley 
172176b2cf59SMatthew Knepley   Input Parameters:
172276b2cf59SMatthew Knepley . snes - The nonlinear solver context
172376b2cf59SMatthew Knepley . func - The function
172476b2cf59SMatthew Knepley 
172576b2cf59SMatthew Knepley   Calling sequence of func:
1726b5d30489SBarry Smith . func (SNES snes, PetscInt step);
172776b2cf59SMatthew Knepley 
172876b2cf59SMatthew Knepley . step - The current step of the iteration
172976b2cf59SMatthew Knepley 
173076b2cf59SMatthew Knepley   Level: intermediate
173176b2cf59SMatthew Knepley 
173276b2cf59SMatthew Knepley .keywords: SNES, update
1733b5d30489SBarry Smith 
173485385478SLisandro Dalcin .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve()
173576b2cf59SMatthew Knepley @*/
173663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
173776b2cf59SMatthew Knepley {
173876b2cf59SMatthew Knepley   PetscFunctionBegin;
17394482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
1740e7788613SBarry Smith   snes->ops->update = func;
174176b2cf59SMatthew Knepley   PetscFunctionReturn(0);
174276b2cf59SMatthew Knepley }
174376b2cf59SMatthew Knepley 
1744e74ef692SMatthew Knepley #undef __FUNCT__
1745e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
174676b2cf59SMatthew Knepley /*@
174776b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
174876b2cf59SMatthew Knepley 
174976b2cf59SMatthew Knepley   Not collective
175076b2cf59SMatthew Knepley 
175176b2cf59SMatthew Knepley   Input Parameters:
175276b2cf59SMatthew Knepley . snes - The nonlinear solver context
175376b2cf59SMatthew Knepley . step - The current step of the iteration
175476b2cf59SMatthew Knepley 
1755205452f4SMatthew Knepley   Level: intermediate
1756205452f4SMatthew Knepley 
175776b2cf59SMatthew Knepley .keywords: SNES, update
1758a6570f20SBarry Smith .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC()
175976b2cf59SMatthew Knepley @*/
176063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
176176b2cf59SMatthew Knepley {
176276b2cf59SMatthew Knepley   PetscFunctionBegin;
176376b2cf59SMatthew Knepley   PetscFunctionReturn(0);
176476b2cf59SMatthew Knepley }
176576b2cf59SMatthew Knepley 
17664a2ae208SSatish Balay #undef __FUNCT__
17674a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
17689b94acceSBarry Smith /*
17699b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
17709b94acceSBarry Smith    positive parameter delta.
17719b94acceSBarry Smith 
17729b94acceSBarry Smith     Input Parameters:
1773c7afd0dbSLois Curfman McInnes +   snes - the SNES context
17749b94acceSBarry Smith .   y - approximate solution of linear system
17759b94acceSBarry Smith .   fnorm - 2-norm of current function
1776c7afd0dbSLois Curfman McInnes -   delta - trust region size
17779b94acceSBarry Smith 
17789b94acceSBarry Smith     Output Parameters:
1779c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
17809b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
17819b94acceSBarry Smith     region, and exceeds zero otherwise.
1782c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
17839b94acceSBarry Smith 
17849b94acceSBarry Smith     Note:
17854b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
17869b94acceSBarry Smith     is set to be the maximum allowable step size.
17879b94acceSBarry Smith 
17889b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
17899b94acceSBarry Smith */
1790dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
17919b94acceSBarry Smith {
1792064f8208SBarry Smith   PetscReal      nrm;
1793ea709b57SSatish Balay   PetscScalar    cnorm;
1794dfbe8321SBarry Smith   PetscErrorCode ierr;
17953a40ed3dSBarry Smith 
17963a40ed3dSBarry Smith   PetscFunctionBegin;
17974482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17984482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1799c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1800184914b5SBarry Smith 
1801064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1802064f8208SBarry Smith   if (nrm > *delta) {
1803064f8208SBarry Smith      nrm = *delta/nrm;
1804064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1805064f8208SBarry Smith      cnorm = nrm;
18062dcb1b2aSMatthew Knepley      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
18079b94acceSBarry Smith      *ynorm = *delta;
18089b94acceSBarry Smith   } else {
18099b94acceSBarry Smith      *gpnorm = 0.0;
1810064f8208SBarry Smith      *ynorm = nrm;
18119b94acceSBarry Smith   }
18123a40ed3dSBarry Smith   PetscFunctionReturn(0);
18139b94acceSBarry Smith }
18149b94acceSBarry Smith 
18154a2ae208SSatish Balay #undef __FUNCT__
18164a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
18176ce558aeSBarry Smith /*@C
1818f69a0ea3SMatthew Knepley    SNESSolve - Solves a nonlinear system F(x) = b.
1819f69a0ea3SMatthew Knepley    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
18209b94acceSBarry Smith 
1821c7afd0dbSLois Curfman McInnes    Collective on SNES
1822c7afd0dbSLois Curfman McInnes 
1823b2002411SLois Curfman McInnes    Input Parameters:
1824c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1825f69a0ea3SMatthew Knepley .  b - the constant part of the equation, or PETSC_NULL to use zero.
182685385478SLisandro Dalcin -  x - the solution vector.
18279b94acceSBarry Smith 
1828b2002411SLois Curfman McInnes    Notes:
18298ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
18308ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
18318ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
18328ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
18338ddd3da0SLois Curfman McInnes 
183436851e7fSLois Curfman McInnes    Level: beginner
183536851e7fSLois Curfman McInnes 
18369b94acceSBarry Smith .keywords: SNES, nonlinear, solve
18379b94acceSBarry Smith 
183885385478SLisandro Dalcin .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian()
18399b94acceSBarry Smith @*/
1840f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
18419b94acceSBarry Smith {
1842dfbe8321SBarry Smith   PetscErrorCode ierr;
1843f1af5d2fSBarry Smith   PetscTruth     flg;
1844eabae89aSBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
1845eabae89aSBarry Smith   PetscViewer    viewer;
1846052efed2SBarry Smith 
18473a40ed3dSBarry Smith   PetscFunctionBegin;
18484482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1849f69a0ea3SMatthew Knepley   PetscValidHeaderSpecific(x,VEC_COOKIE,3);
1850f69a0ea3SMatthew Knepley   PetscCheckSameComm(snes,1,x,3);
185185385478SLisandro Dalcin   if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2);
185285385478SLisandro Dalcin   if (b) PetscCheckSameComm(snes,1,b,2);
185385385478SLisandro Dalcin 
185485385478SLisandro Dalcin   /* set solution vector */
185585385478SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr);
185685385478SLisandro Dalcin   if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); }
185785385478SLisandro Dalcin   snes->vec_sol = x;
185885385478SLisandro Dalcin   /* set afine vector if provided */
185985385478SLisandro Dalcin   if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); }
186085385478SLisandro Dalcin   if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); }
186185385478SLisandro Dalcin   snes->vec_rhs = b;
186285385478SLisandro Dalcin 
186385385478SLisandro Dalcin   if (!snes->vec_func && snes->vec_rhs) {
186485385478SLisandro Dalcin     ierr = VecDuplicate(b, &snes->vec_func); CHKERRQ(ierr);
186570e92668SMatthew Knepley   }
18663f149594SLisandro Dalcin 
186770e92668SMatthew Knepley   ierr = SNESSetUp(snes);CHKERRQ(ierr);
18683f149594SLisandro Dalcin 
1869abc0a331SBarry Smith   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
187050ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1871d5e45103SBarry Smith 
18723f149594SLisandro Dalcin   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
18733f149594SLisandro Dalcin 
1874e7788613SBarry Smith   ierr = PetscExceptionTry1((*(snes)->ops->solve)(snes),PETSC_ERR_ARG_DOMAIN);
1875d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
187638f152feSBarry Smith     /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */
187719717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr);
1878c671dbe7SSatish Balay   } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) {
1879d5e45103SBarry Smith     /* translate exception into SNES not converged reason */
1880d5e45103SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
188138f152feSBarry Smith     ierr = 0;
188238f152feSBarry Smith   }
188338f152feSBarry Smith   CHKERRQ(ierr);
188485385478SLisandro Dalcin   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
188585385478SLisandro Dalcin 
18863f149594SLisandro Dalcin   if (!snes->reason) {
18873f149594SLisandro Dalcin     SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
18883f149594SLisandro Dalcin   }
18893f149594SLisandro Dalcin 
18907adad957SLisandro Dalcin   ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1891eabae89aSBarry Smith   if (flg && !PetscPreLoadingOn) {
18927adad957SLisandro Dalcin     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr);
1893eabae89aSBarry Smith     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
1894eabae89aSBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
1895eabae89aSBarry Smith   }
1896eabae89aSBarry Smith 
18977adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1898da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
18995968eb51SBarry Smith   if (snes->printreason) {
19005968eb51SBarry Smith     if (snes->reason > 0) {
19017adad957SLisandro Dalcin       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
19025968eb51SBarry Smith     } else {
19037adad957SLisandro Dalcin       ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
19045968eb51SBarry Smith     }
19055968eb51SBarry Smith   }
19065968eb51SBarry Smith 
19073a40ed3dSBarry Smith   PetscFunctionReturn(0);
19089b94acceSBarry Smith }
19099b94acceSBarry Smith 
19109b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
19119b94acceSBarry Smith 
19124a2ae208SSatish Balay #undef __FUNCT__
19134a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
191482bf6240SBarry Smith /*@C
19154b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
19169b94acceSBarry Smith 
1917fee21e36SBarry Smith    Collective on SNES
1918fee21e36SBarry Smith 
1919c7afd0dbSLois Curfman McInnes    Input Parameters:
1920c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1921454a90a3SBarry Smith -  type - a known method
1922c7afd0dbSLois Curfman McInnes 
1923c7afd0dbSLois Curfman McInnes    Options Database Key:
1924454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1925c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1926ae12b187SLois Curfman McInnes 
19279b94acceSBarry Smith    Notes:
1928e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
19294b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
1930c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
19314b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
1932c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
19339b94acceSBarry Smith 
1934ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1935ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1936ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1937ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1938ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1939ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1940ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1941ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1942ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1943b0a32e0cSBarry Smith   appropriate method.
194436851e7fSLois Curfman McInnes 
194536851e7fSLois Curfman McInnes   Level: intermediate
1946a703fe33SLois Curfman McInnes 
1947454a90a3SBarry Smith .keywords: SNES, set, type
1948435da068SBarry Smith 
1949435da068SBarry Smith .seealso: SNESType, SNESCreate()
1950435da068SBarry Smith 
19519b94acceSBarry Smith @*/
1952e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type)
19539b94acceSBarry Smith {
1954dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
19556831982aSBarry Smith   PetscTruth     match;
19563a40ed3dSBarry Smith 
19573a40ed3dSBarry Smith   PetscFunctionBegin;
19584482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19594482741eSBarry Smith   PetscValidCharPointer(type,2);
196082bf6240SBarry Smith 
19616831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
19620f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
196392ff6ae8SBarry Smith 
19647adad957SLisandro Dalcin   ierr =  PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr);
1965958c9bccSBarry Smith   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
196675396ef9SLisandro Dalcin   /* Destroy the previous private SNES context */
196775396ef9SLisandro Dalcin   if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); }
196875396ef9SLisandro Dalcin   /* Reinitialize function pointers in SNESOps structure */
196975396ef9SLisandro Dalcin   snes->ops->setup          = 0;
197075396ef9SLisandro Dalcin   snes->ops->solve          = 0;
197175396ef9SLisandro Dalcin   snes->ops->view           = 0;
197275396ef9SLisandro Dalcin   snes->ops->setfromoptions = 0;
197375396ef9SLisandro Dalcin   snes->ops->destroy        = 0;
197475396ef9SLisandro Dalcin   /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */
197575396ef9SLisandro Dalcin   snes->setupcalled = PETSC_FALSE;
19763a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1977454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
19783a40ed3dSBarry Smith   PetscFunctionReturn(0);
19799b94acceSBarry Smith }
19809b94acceSBarry Smith 
1981a847f771SSatish Balay 
19829b94acceSBarry Smith /* --------------------------------------------------------------------- */
19834a2ae208SSatish Balay #undef __FUNCT__
19844a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
198552baeb72SSatish Balay /*@
19869b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1987f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
19889b94acceSBarry Smith 
1989fee21e36SBarry Smith    Not Collective
1990fee21e36SBarry Smith 
199136851e7fSLois Curfman McInnes    Level: advanced
199236851e7fSLois Curfman McInnes 
19939b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
19949b94acceSBarry Smith 
19959b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
19969b94acceSBarry Smith @*/
199763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
19989b94acceSBarry Smith {
1999dfbe8321SBarry Smith   PetscErrorCode ierr;
200082bf6240SBarry Smith 
20013a40ed3dSBarry Smith   PetscFunctionBegin;
20021441b1d3SBarry Smith   ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
20034c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
20043a40ed3dSBarry Smith   PetscFunctionReturn(0);
20059b94acceSBarry Smith }
20069b94acceSBarry Smith 
20074a2ae208SSatish Balay #undef __FUNCT__
20084a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
20099b94acceSBarry Smith /*@C
20109a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
20119b94acceSBarry Smith 
2012c7afd0dbSLois Curfman McInnes    Not Collective
2013c7afd0dbSLois Curfman McInnes 
20149b94acceSBarry Smith    Input Parameter:
20154b0e389bSBarry Smith .  snes - nonlinear solver context
20169b94acceSBarry Smith 
20179b94acceSBarry Smith    Output Parameter:
20183a7fca6bSBarry Smith .  type - SNES method (a character string)
20199b94acceSBarry Smith 
202036851e7fSLois Curfman McInnes    Level: intermediate
202136851e7fSLois Curfman McInnes 
2022454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
20239b94acceSBarry Smith @*/
202463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type)
20259b94acceSBarry Smith {
20263a40ed3dSBarry Smith   PetscFunctionBegin;
20274482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
20284482741eSBarry Smith   PetscValidPointer(type,2);
20297adad957SLisandro Dalcin   *type = ((PetscObject)snes)->type_name;
20303a40ed3dSBarry Smith   PetscFunctionReturn(0);
20319b94acceSBarry Smith }
20329b94acceSBarry Smith 
20334a2ae208SSatish Balay #undef __FUNCT__
20344a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
203552baeb72SSatish Balay /*@
20369b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
20379b94acceSBarry Smith    stored.
20389b94acceSBarry Smith 
2039c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2040c7afd0dbSLois Curfman McInnes 
20419b94acceSBarry Smith    Input Parameter:
20429b94acceSBarry Smith .  snes - the SNES context
20439b94acceSBarry Smith 
20449b94acceSBarry Smith    Output Parameter:
20459b94acceSBarry Smith .  x - the solution
20469b94acceSBarry Smith 
204770e92668SMatthew Knepley    Level: intermediate
204836851e7fSLois Curfman McInnes 
20499b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
20509b94acceSBarry Smith 
205185385478SLisandro Dalcin .seealso:  SNESGetSolutionUpdate(), SNESGetFunction()
20529b94acceSBarry Smith @*/
205363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
20549b94acceSBarry Smith {
20553a40ed3dSBarry Smith   PetscFunctionBegin;
20564482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
20574482741eSBarry Smith   PetscValidPointer(x,2);
205885385478SLisandro Dalcin   *x = snes->vec_sol;
205970e92668SMatthew Knepley   PetscFunctionReturn(0);
206070e92668SMatthew Knepley }
206170e92668SMatthew Knepley 
206270e92668SMatthew Knepley #undef __FUNCT__
20634a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
206452baeb72SSatish Balay /*@
20659b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
20669b94acceSBarry Smith    stored.
20679b94acceSBarry Smith 
2068c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2069c7afd0dbSLois Curfman McInnes 
20709b94acceSBarry Smith    Input Parameter:
20719b94acceSBarry Smith .  snes - the SNES context
20729b94acceSBarry Smith 
20739b94acceSBarry Smith    Output Parameter:
20749b94acceSBarry Smith .  x - the solution update
20759b94acceSBarry Smith 
207636851e7fSLois Curfman McInnes    Level: advanced
207736851e7fSLois Curfman McInnes 
20789b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
20799b94acceSBarry Smith 
208085385478SLisandro Dalcin .seealso: SNESGetSolution(), SNESGetFunction()
20819b94acceSBarry Smith @*/
208263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
20839b94acceSBarry Smith {
20843a40ed3dSBarry Smith   PetscFunctionBegin;
20854482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
20864482741eSBarry Smith   PetscValidPointer(x,2);
208785385478SLisandro Dalcin   *x = snes->vec_sol_update;
20883a40ed3dSBarry Smith   PetscFunctionReturn(0);
20899b94acceSBarry Smith }
20909b94acceSBarry Smith 
20914a2ae208SSatish Balay #undef __FUNCT__
20924a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
20939b94acceSBarry Smith /*@C
20943638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
20959b94acceSBarry Smith 
2096c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2097c7afd0dbSLois Curfman McInnes 
20989b94acceSBarry Smith    Input Parameter:
20999b94acceSBarry Smith .  snes - the SNES context
21009b94acceSBarry Smith 
21019b94acceSBarry Smith    Output Parameter:
21027bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
210370e92668SMatthew Knepley .  func - the function (or PETSC_NULL)
210470e92668SMatthew Knepley -  ctx - the function context (or PETSC_NULL)
21059b94acceSBarry Smith 
210636851e7fSLois Curfman McInnes    Level: advanced
210736851e7fSLois Curfman McInnes 
2108a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
21099b94acceSBarry Smith 
21104b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
21119b94acceSBarry Smith @*/
211270e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
21139b94acceSBarry Smith {
21143a40ed3dSBarry Smith   PetscFunctionBegin;
21154482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
211685385478SLisandro Dalcin   if (r)    *r    = snes->vec_func;
2117e7788613SBarry Smith   if (func) *func = snes->ops->computefunction;
211870e92668SMatthew Knepley   if (ctx)  *ctx  = snes->funP;
21193a40ed3dSBarry Smith   PetscFunctionReturn(0);
21209b94acceSBarry Smith }
21219b94acceSBarry Smith 
21224a2ae208SSatish Balay #undef __FUNCT__
21234a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
21243c7409f5SSatish Balay /*@C
21253c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2126d850072dSLois Curfman McInnes    SNES options in the database.
21273c7409f5SSatish Balay 
2128fee21e36SBarry Smith    Collective on SNES
2129fee21e36SBarry Smith 
2130c7afd0dbSLois Curfman McInnes    Input Parameter:
2131c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2132c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2133c7afd0dbSLois Curfman McInnes 
2134d850072dSLois Curfman McInnes    Notes:
2135a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2136c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2137d850072dSLois Curfman McInnes 
213836851e7fSLois Curfman McInnes    Level: advanced
213936851e7fSLois Curfman McInnes 
21403c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2141a86d99e1SLois Curfman McInnes 
2142a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
21433c7409f5SSatish Balay @*/
214463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
21453c7409f5SSatish Balay {
2146dfbe8321SBarry Smith   PetscErrorCode ierr;
21473c7409f5SSatish Balay 
21483a40ed3dSBarry Smith   PetscFunctionBegin;
21494482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2150639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
215194b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
21523a40ed3dSBarry Smith   PetscFunctionReturn(0);
21533c7409f5SSatish Balay }
21543c7409f5SSatish Balay 
21554a2ae208SSatish Balay #undef __FUNCT__
21564a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
21573c7409f5SSatish Balay /*@C
2158f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2159d850072dSLois Curfman McInnes    SNES options in the database.
21603c7409f5SSatish Balay 
2161fee21e36SBarry Smith    Collective on SNES
2162fee21e36SBarry Smith 
2163c7afd0dbSLois Curfman McInnes    Input Parameters:
2164c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2165c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2166c7afd0dbSLois Curfman McInnes 
2167d850072dSLois Curfman McInnes    Notes:
2168a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2169c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2170d850072dSLois Curfman McInnes 
217136851e7fSLois Curfman McInnes    Level: advanced
217236851e7fSLois Curfman McInnes 
21733c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2174a86d99e1SLois Curfman McInnes 
2175a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
21763c7409f5SSatish Balay @*/
217763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
21783c7409f5SSatish Balay {
2179dfbe8321SBarry Smith   PetscErrorCode ierr;
21803c7409f5SSatish Balay 
21813a40ed3dSBarry Smith   PetscFunctionBegin;
21824482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2183639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
218494b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
21853a40ed3dSBarry Smith   PetscFunctionReturn(0);
21863c7409f5SSatish Balay }
21873c7409f5SSatish Balay 
21884a2ae208SSatish Balay #undef __FUNCT__
21894a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
21909ab63eb5SSatish Balay /*@C
21913c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
21923c7409f5SSatish Balay    SNES options in the database.
21933c7409f5SSatish Balay 
2194c7afd0dbSLois Curfman McInnes    Not Collective
2195c7afd0dbSLois Curfman McInnes 
21963c7409f5SSatish Balay    Input Parameter:
21973c7409f5SSatish Balay .  snes - the SNES context
21983c7409f5SSatish Balay 
21993c7409f5SSatish Balay    Output Parameter:
22003c7409f5SSatish Balay .  prefix - pointer to the prefix string used
22013c7409f5SSatish Balay 
22029ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
22039ab63eb5SSatish Balay    sufficient length to hold the prefix.
22049ab63eb5SSatish Balay 
220536851e7fSLois Curfman McInnes    Level: advanced
220636851e7fSLois Curfman McInnes 
22073c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2208a86d99e1SLois Curfman McInnes 
2209a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
22103c7409f5SSatish Balay @*/
2211e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
22123c7409f5SSatish Balay {
2213dfbe8321SBarry Smith   PetscErrorCode ierr;
22143c7409f5SSatish Balay 
22153a40ed3dSBarry Smith   PetscFunctionBegin;
22164482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2217639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
22183a40ed3dSBarry Smith   PetscFunctionReturn(0);
22193c7409f5SSatish Balay }
22203c7409f5SSatish Balay 
2221b2002411SLois Curfman McInnes 
22224a2ae208SSatish Balay #undef __FUNCT__
22234a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
22243cea93caSBarry Smith /*@C
22253cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
22263cea93caSBarry Smith 
22277f6c08e0SMatthew Knepley   Level: advanced
22283cea93caSBarry Smith @*/
222963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2230b2002411SLois Curfman McInnes {
2231e2d1d2b7SBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
2232dfbe8321SBarry Smith   PetscErrorCode ierr;
2233b2002411SLois Curfman McInnes 
2234b2002411SLois Curfman McInnes   PetscFunctionBegin;
2235b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2236c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2237b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2238b2002411SLois Curfman McInnes }
2239da9b6338SBarry Smith 
2240da9b6338SBarry Smith #undef __FUNCT__
2241da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
224263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2243da9b6338SBarry Smith {
2244dfbe8321SBarry Smith   PetscErrorCode ierr;
224577431f27SBarry Smith   PetscInt       N,i,j;
2246da9b6338SBarry Smith   Vec            u,uh,fh;
2247da9b6338SBarry Smith   PetscScalar    value;
2248da9b6338SBarry Smith   PetscReal      norm;
2249da9b6338SBarry Smith 
2250da9b6338SBarry Smith   PetscFunctionBegin;
2251da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2252da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2253da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2254da9b6338SBarry Smith 
2255da9b6338SBarry Smith   /* currently only works for sequential */
2256da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2257da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2258da9b6338SBarry Smith   for (i=0; i<N; i++) {
2259da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
226077431f27SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2261da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2262ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2263da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
22643ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2265da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
226677431f27SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2267da9b6338SBarry Smith       value = -value;
2268da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2269da9b6338SBarry Smith     }
2270da9b6338SBarry Smith   }
2271da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2272da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2273da9b6338SBarry Smith   PetscFunctionReturn(0);
2274da9b6338SBarry Smith }
227571f87433Sdalcinl 
227671f87433Sdalcinl #undef __FUNCT__
2277fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPSetUseEW"
227871f87433Sdalcinl /*@
2279fa9f3622SBarry Smith    SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
228071f87433Sdalcinl    computing relative tolerance for linear solvers within an inexact
228171f87433Sdalcinl    Newton method.
228271f87433Sdalcinl 
228371f87433Sdalcinl    Collective on SNES
228471f87433Sdalcinl 
228571f87433Sdalcinl    Input Parameters:
228671f87433Sdalcinl +  snes - SNES context
228771f87433Sdalcinl -  flag - PETSC_TRUE or PETSC_FALSE
228871f87433Sdalcinl 
228971f87433Sdalcinl    Notes:
229071f87433Sdalcinl    Currently, the default is to use a constant relative tolerance for
229171f87433Sdalcinl    the inner linear solvers.  Alternatively, one can use the
229271f87433Sdalcinl    Eisenstat-Walker method, where the relative convergence tolerance
229371f87433Sdalcinl    is reset at each Newton iteration according progress of the nonlinear
229471f87433Sdalcinl    solver.
229571f87433Sdalcinl 
229671f87433Sdalcinl    Level: advanced
229771f87433Sdalcinl 
229871f87433Sdalcinl    Reference:
229971f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
230071f87433Sdalcinl    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
230171f87433Sdalcinl 
230271f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
230371f87433Sdalcinl 
2304fa9f3622SBarry Smith .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
230571f87433Sdalcinl @*/
2306fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag)
230771f87433Sdalcinl {
230871f87433Sdalcinl   PetscFunctionBegin;
230971f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
231071f87433Sdalcinl   snes->ksp_ewconv = flag;
231171f87433Sdalcinl   PetscFunctionReturn(0);
231271f87433Sdalcinl }
231371f87433Sdalcinl 
231471f87433Sdalcinl #undef __FUNCT__
2315fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPGetUseEW"
231671f87433Sdalcinl /*@
2317fa9f3622SBarry Smith    SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method
231871f87433Sdalcinl    for computing relative tolerance for linear solvers within an
231971f87433Sdalcinl    inexact Newton method.
232071f87433Sdalcinl 
232171f87433Sdalcinl    Not Collective
232271f87433Sdalcinl 
232371f87433Sdalcinl    Input Parameter:
232471f87433Sdalcinl .  snes - SNES context
232571f87433Sdalcinl 
232671f87433Sdalcinl    Output Parameter:
232771f87433Sdalcinl .  flag - PETSC_TRUE or PETSC_FALSE
232871f87433Sdalcinl 
232971f87433Sdalcinl    Notes:
233071f87433Sdalcinl    Currently, the default is to use a constant relative tolerance for
233171f87433Sdalcinl    the inner linear solvers.  Alternatively, one can use the
233271f87433Sdalcinl    Eisenstat-Walker method, where the relative convergence tolerance
233371f87433Sdalcinl    is reset at each Newton iteration according progress of the nonlinear
233471f87433Sdalcinl    solver.
233571f87433Sdalcinl 
233671f87433Sdalcinl    Level: advanced
233771f87433Sdalcinl 
233871f87433Sdalcinl    Reference:
233971f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
234071f87433Sdalcinl    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
234171f87433Sdalcinl 
234271f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton
234371f87433Sdalcinl 
2344fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW()
234571f87433Sdalcinl @*/
2346fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag)
234771f87433Sdalcinl {
234871f87433Sdalcinl   PetscFunctionBegin;
234971f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
235071f87433Sdalcinl   PetscValidPointer(flag,2);
235171f87433Sdalcinl   *flag = snes->ksp_ewconv;
235271f87433Sdalcinl   PetscFunctionReturn(0);
235371f87433Sdalcinl }
235471f87433Sdalcinl 
235571f87433Sdalcinl #undef __FUNCT__
2356fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPSetParametersEW"
235771f87433Sdalcinl /*@
2358fa9f3622SBarry Smith    SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker
235971f87433Sdalcinl    convergence criteria for the linear solvers within an inexact
236071f87433Sdalcinl    Newton method.
236171f87433Sdalcinl 
236271f87433Sdalcinl    Collective on SNES
236371f87433Sdalcinl 
236471f87433Sdalcinl    Input Parameters:
236571f87433Sdalcinl +    snes - SNES context
236671f87433Sdalcinl .    version - version 1, 2 (default is 2) or 3
236771f87433Sdalcinl .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
236871f87433Sdalcinl .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
236971f87433Sdalcinl .    gamma - multiplicative factor for version 2 rtol computation
237071f87433Sdalcinl              (0 <= gamma2 <= 1)
237171f87433Sdalcinl .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
237271f87433Sdalcinl .    alpha2 - power for safeguard
237371f87433Sdalcinl -    threshold - threshold for imposing safeguard (0 < threshold < 1)
237471f87433Sdalcinl 
237571f87433Sdalcinl    Note:
237671f87433Sdalcinl    Version 3 was contributed by Luis Chacon, June 2006.
237771f87433Sdalcinl 
237871f87433Sdalcinl    Use PETSC_DEFAULT to retain the default for any of the parameters.
237971f87433Sdalcinl 
238071f87433Sdalcinl    Level: advanced
238171f87433Sdalcinl 
238271f87433Sdalcinl    Reference:
238371f87433Sdalcinl    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
238471f87433Sdalcinl    inexact Newton method", Utah State University Math. Stat. Dept. Res.
238571f87433Sdalcinl    Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput.
238671f87433Sdalcinl 
238771f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, set, parameters
238871f87433Sdalcinl 
2389fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW()
239071f87433Sdalcinl @*/
2391fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max,
239271f87433Sdalcinl 							    PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold)
239371f87433Sdalcinl {
2394fa9f3622SBarry Smith   SNESKSPEW *kctx;
239571f87433Sdalcinl   PetscFunctionBegin;
239671f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2397fa9f3622SBarry Smith   kctx = (SNESKSPEW*)snes->kspconvctx;
239871f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
239971f87433Sdalcinl 
240071f87433Sdalcinl   if (version != PETSC_DEFAULT)   kctx->version   = version;
240171f87433Sdalcinl   if (rtol_0 != PETSC_DEFAULT)    kctx->rtol_0    = rtol_0;
240271f87433Sdalcinl   if (rtol_max != PETSC_DEFAULT)  kctx->rtol_max  = rtol_max;
240371f87433Sdalcinl   if (gamma != PETSC_DEFAULT)     kctx->gamma     = gamma;
240471f87433Sdalcinl   if (alpha != PETSC_DEFAULT)     kctx->alpha     = alpha;
240571f87433Sdalcinl   if (alpha2 != PETSC_DEFAULT)    kctx->alpha2    = alpha2;
240671f87433Sdalcinl   if (threshold != PETSC_DEFAULT) kctx->threshold = threshold;
240771f87433Sdalcinl 
240871f87433Sdalcinl   if (kctx->version < 1 || kctx->version > 3) {
240971f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version);
241071f87433Sdalcinl   }
241171f87433Sdalcinl   if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) {
241271f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0);
241371f87433Sdalcinl   }
241471f87433Sdalcinl   if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) {
241571f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max);
241671f87433Sdalcinl   }
241771f87433Sdalcinl   if (kctx->gamma < 0.0 || kctx->gamma > 1.0) {
241871f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma);
241971f87433Sdalcinl   }
242071f87433Sdalcinl   if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) {
242171f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha);
242271f87433Sdalcinl   }
242371f87433Sdalcinl   if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) {
242471f87433Sdalcinl     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold);
242571f87433Sdalcinl   }
242671f87433Sdalcinl   PetscFunctionReturn(0);
242771f87433Sdalcinl }
242871f87433Sdalcinl 
242971f87433Sdalcinl #undef __FUNCT__
2430fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPGetParametersEW"
243171f87433Sdalcinl /*@
2432fa9f3622SBarry Smith    SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker
243371f87433Sdalcinl    convergence criteria for the linear solvers within an inexact
243471f87433Sdalcinl    Newton method.
243571f87433Sdalcinl 
243671f87433Sdalcinl    Not Collective
243771f87433Sdalcinl 
243871f87433Sdalcinl    Input Parameters:
243971f87433Sdalcinl      snes - SNES context
244071f87433Sdalcinl 
244171f87433Sdalcinl    Output Parameters:
244271f87433Sdalcinl +    version - version 1, 2 (default is 2) or 3
244371f87433Sdalcinl .    rtol_0 - initial relative tolerance (0 <= rtol_0 < 1)
244471f87433Sdalcinl .    rtol_max - maximum relative tolerance (0 <= rtol_max < 1)
244571f87433Sdalcinl .    gamma - multiplicative factor for version 2 rtol computation
244671f87433Sdalcinl              (0 <= gamma2 <= 1)
244771f87433Sdalcinl .    alpha - power for version 2 rtol computation (1 < alpha <= 2)
244871f87433Sdalcinl .    alpha2 - power for safeguard
244971f87433Sdalcinl -    threshold - threshold for imposing safeguard (0 < threshold < 1)
245071f87433Sdalcinl 
245171f87433Sdalcinl    Level: advanced
245271f87433Sdalcinl 
245371f87433Sdalcinl .keywords: SNES, KSP, Eisenstat, Walker, get, parameters
245471f87433Sdalcinl 
2455fa9f3622SBarry Smith .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW()
245671f87433Sdalcinl @*/
2457fa9f3622SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max,
245871f87433Sdalcinl 							    PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold)
245971f87433Sdalcinl {
2460fa9f3622SBarry Smith   SNESKSPEW *kctx;
246171f87433Sdalcinl   PetscFunctionBegin;
246271f87433Sdalcinl   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2463fa9f3622SBarry Smith   kctx = (SNESKSPEW*)snes->kspconvctx;
246471f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing");
246571f87433Sdalcinl   if(version)   *version   = kctx->version;
246671f87433Sdalcinl   if(rtol_0)    *rtol_0    = kctx->rtol_0;
246771f87433Sdalcinl   if(rtol_max)  *rtol_max  = kctx->rtol_max;
246871f87433Sdalcinl   if(gamma)     *gamma     = kctx->gamma;
246971f87433Sdalcinl   if(alpha)     *alpha     = kctx->alpha;
247071f87433Sdalcinl   if(alpha2)    *alpha2    = kctx->alpha2;
247171f87433Sdalcinl   if(threshold) *threshold = kctx->threshold;
247271f87433Sdalcinl   PetscFunctionReturn(0);
247371f87433Sdalcinl }
247471f87433Sdalcinl 
247571f87433Sdalcinl #undef __FUNCT__
2476fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPEW_PreSolve"
2477fa9f3622SBarry Smith static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x)
247871f87433Sdalcinl {
247971f87433Sdalcinl   PetscErrorCode ierr;
2480fa9f3622SBarry Smith   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
248171f87433Sdalcinl   PetscReal      rtol=PETSC_DEFAULT,stol;
248271f87433Sdalcinl 
248371f87433Sdalcinl   PetscFunctionBegin;
248471f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
248571f87433Sdalcinl   if (!snes->iter) { /* first time in, so use the original user rtol */
248671f87433Sdalcinl     rtol = kctx->rtol_0;
248771f87433Sdalcinl   } else {
248871f87433Sdalcinl     if (kctx->version == 1) {
248971f87433Sdalcinl       rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last;
249071f87433Sdalcinl       if (rtol < 0.0) rtol = -rtol;
249171f87433Sdalcinl       stol = pow(kctx->rtol_last,kctx->alpha2);
249271f87433Sdalcinl       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
249371f87433Sdalcinl     } else if (kctx->version == 2) {
249471f87433Sdalcinl       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
249571f87433Sdalcinl       stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha);
249671f87433Sdalcinl       if (stol > kctx->threshold) rtol = PetscMax(rtol,stol);
249771f87433Sdalcinl     } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */
249871f87433Sdalcinl       rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha);
249971f87433Sdalcinl       /* safeguard: avoid sharp decrease of rtol */
250071f87433Sdalcinl       stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha);
250171f87433Sdalcinl       stol = PetscMax(rtol,stol);
250271f87433Sdalcinl       rtol = PetscMin(kctx->rtol_0,stol);
250371f87433Sdalcinl       /* safeguard: avoid oversolving */
250471f87433Sdalcinl       stol = kctx->gamma*(snes->ttol)/snes->norm;
250571f87433Sdalcinl       stol = PetscMax(rtol,stol);
250671f87433Sdalcinl       rtol = PetscMin(kctx->rtol_0,stol);
250771f87433Sdalcinl     } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version);
250871f87433Sdalcinl   }
250971f87433Sdalcinl   /* safeguard: avoid rtol greater than one */
251071f87433Sdalcinl   rtol = PetscMin(rtol,kctx->rtol_max);
251171f87433Sdalcinl   ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
251271f87433Sdalcinl   ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr);
251371f87433Sdalcinl   PetscFunctionReturn(0);
251471f87433Sdalcinl }
251571f87433Sdalcinl 
251671f87433Sdalcinl #undef __FUNCT__
2517fa9f3622SBarry Smith #define __FUNCT__ "SNESKSPEW_PostSolve"
2518fa9f3622SBarry Smith static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x)
251971f87433Sdalcinl {
252071f87433Sdalcinl   PetscErrorCode ierr;
2521fa9f3622SBarry Smith   SNESKSPEW      *kctx = (SNESKSPEW*)snes->kspconvctx;
252271f87433Sdalcinl   PCSide         pcside;
252371f87433Sdalcinl   Vec            lres;
252471f87433Sdalcinl 
252571f87433Sdalcinl   PetscFunctionBegin;
252671f87433Sdalcinl   if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists");
252771f87433Sdalcinl   ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr);
252871f87433Sdalcinl   ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr);
252971f87433Sdalcinl   if (kctx->version == 1) {
253071f87433Sdalcinl     ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr);
253171f87433Sdalcinl     if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */
253271f87433Sdalcinl       /* KSP residual is true linear residual */
253371f87433Sdalcinl       ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr);
253471f87433Sdalcinl     } else {
253571f87433Sdalcinl       /* KSP residual is preconditioned residual */
253671f87433Sdalcinl       /* compute true linear residual norm */
253771f87433Sdalcinl       ierr = VecDuplicate(b,&lres);CHKERRQ(ierr);
253871f87433Sdalcinl       ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr);
253971f87433Sdalcinl       ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr);
254071f87433Sdalcinl       ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr);
254171f87433Sdalcinl       ierr = VecDestroy(lres);CHKERRQ(ierr);
254271f87433Sdalcinl     }
254371f87433Sdalcinl   }
254471f87433Sdalcinl   PetscFunctionReturn(0);
254571f87433Sdalcinl }
254671f87433Sdalcinl 
254771f87433Sdalcinl #undef __FUNCT__
254871f87433Sdalcinl #define __FUNCT__ "SNES_KSPSolve"
254971f87433Sdalcinl PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x)
255071f87433Sdalcinl {
255171f87433Sdalcinl   PetscErrorCode ierr;
255271f87433Sdalcinl 
255371f87433Sdalcinl   PetscFunctionBegin;
2554fa9f3622SBarry Smith   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr);  }
255571f87433Sdalcinl   ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);
2556fa9f3622SBarry Smith   if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); }
255771f87433Sdalcinl   PetscFunctionReturn(0);
255871f87433Sdalcinl }
2559