xref: /petsc/src/snes/interface/snes.c (revision ae15b995b5732fffd2de5a75cf61ef7190c6fef1)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
29b94acceSBarry Smith 
3e090d566SSatish Balay #include "src/snes/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 {
459b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *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);
53b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);
544482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
55c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,viewer,2);
5674679c65SBarry Smith 
5732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
58b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
5932077d6dSBarry Smith   if (iascii) {
603a7fca6bSBarry Smith     if (snes->prefix) {
613a7fca6bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr);
623a7fca6bSBarry Smith     } else {
63b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
643a7fca6bSBarry Smith     }
65454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
66454a90a3SBarry Smith     if (type) {
67b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
68184914b5SBarry Smith     } else {
69b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
70184914b5SBarry Smith     }
710ef38995SBarry Smith     if (snes->view) {
72b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
730ef38995SBarry Smith       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
74b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
750ef38995SBarry Smith     }
7677431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
77a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%G, absolute=%G, solution=%G\n",
7870441072SBarry Smith                  snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr);
7977431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr);
8077431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr);
819b94acceSBarry Smith     if (snes->ksp_ewconv) {
829b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
839b94acceSBarry Smith       if (kctx) {
8477431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr);
85a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
86a83599f4SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
879b94acceSBarry Smith       }
889b94acceSBarry Smith     }
890f5bd95cSBarry Smith   } else if (isstring) {
90454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
91b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
9219bcc07fSBarry Smith   }
9394b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
94b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
9594b7f48cSBarry Smith   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
96b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
973a40ed3dSBarry Smith   PetscFunctionReturn(0);
989b94acceSBarry Smith }
999b94acceSBarry Smith 
10076b2cf59SMatthew Knepley /*
10176b2cf59SMatthew Knepley   We retain a list of functions that also take SNES command
10276b2cf59SMatthew Knepley   line options. These are called at the end SNESSetFromOptions()
10376b2cf59SMatthew Knepley */
10476b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5
105a7cc72afSBarry Smith static PetscInt numberofsetfromoptions;
1066849ba73SBarry Smith static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
10776b2cf59SMatthew Knepley 
108e74ef692SMatthew Knepley #undef __FUNCT__
109e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker"
110ac226902SBarry Smith /*@C
11176b2cf59SMatthew Knepley   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
11276b2cf59SMatthew Knepley 
11376b2cf59SMatthew Knepley   Not Collective
11476b2cf59SMatthew Knepley 
11576b2cf59SMatthew Knepley   Input Parameter:
11676b2cf59SMatthew Knepley . snescheck - function that checks for options
11776b2cf59SMatthew Knepley 
11876b2cf59SMatthew Knepley   Level: developer
11976b2cf59SMatthew Knepley 
12076b2cf59SMatthew Knepley .seealso: SNESSetFromOptions()
12176b2cf59SMatthew Knepley @*/
12263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES))
12376b2cf59SMatthew Knepley {
12476b2cf59SMatthew Knepley   PetscFunctionBegin;
12576b2cf59SMatthew Knepley   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
12677431f27SBarry Smith     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS);
12776b2cf59SMatthew Knepley   }
12876b2cf59SMatthew Knepley   othersetfromoptions[numberofsetfromoptions++] = snescheck;
12976b2cf59SMatthew Knepley   PetscFunctionReturn(0);
13076b2cf59SMatthew Knepley }
13176b2cf59SMatthew Knepley 
1324a2ae208SSatish Balay #undef __FUNCT__
1334a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions"
1349b94acceSBarry Smith /*@
13594b7f48cSBarry Smith    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
1369b94acceSBarry Smith 
137c7afd0dbSLois Curfman McInnes    Collective on SNES
138c7afd0dbSLois Curfman McInnes 
1399b94acceSBarry Smith    Input Parameter:
1409b94acceSBarry Smith .  snes - the SNES context
1419b94acceSBarry Smith 
14236851e7fSLois Curfman McInnes    Options Database Keys:
1436831982aSBarry Smith +  -snes_type <type> - ls, tr, umls, umtr, test
14482738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
14582738288SBarry Smith                 of the change in the solution between steps
14670441072SBarry Smith .  -snes_atol <abstol> - absolute tolerance of residual norm
147b39c3a46SLois Curfman McInnes .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
148b39c3a46SLois Curfman McInnes .  -snes_max_it <max_it> - maximum number of iterations
149b39c3a46SLois Curfman McInnes .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
15050ffb88aSMatthew Knepley .  -snes_max_fail <max_fail> - maximum number of failures
151b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
1522492ecdbSBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear
15382738288SBarry Smith                                solver; hence iterations will continue until max_it
1541fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
15582738288SBarry Smith                                of convergence test
156e8105e01SRichard Katz .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
157e8105e01SRichard Katz                                        filename given prints to stdout
158d132466eSBarry Smith .  -snes_vecmonitor - plots solution at each iteration
159d132466eSBarry Smith .  -snes_vecmonitor_update - plots update to solution at each iteration
16082738288SBarry Smith .  -snes_xmonitor - plots residual norm at each iteration
161e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
1625968eb51SBarry Smith .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
1635968eb51SBarry Smith -  -snes_print_converged_reason - print the reason for convergence/divergence after each solve
16482738288SBarry Smith 
16582738288SBarry Smith     Options Database for Eisenstat-Walker method:
1664b27c08aSLois Curfman McInnes +  -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence
1674b27c08aSLois Curfman McInnes .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
16836851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
16936851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
17036851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
17136851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
17236851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
17336851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
17482738288SBarry Smith 
17511ca99fdSLois Curfman McInnes    Notes:
17611ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
17711ca99fdSLois Curfman McInnes    the users manual.
17883e2fdc7SBarry Smith 
17936851e7fSLois Curfman McInnes    Level: beginner
18036851e7fSLois Curfman McInnes 
1819b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1829b94acceSBarry Smith 
18369ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
1849b94acceSBarry Smith @*/
18563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes)
1869b94acceSBarry Smith {
18794b7f48cSBarry Smith   KSP                 ksp;
188186905e3SBarry Smith   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
189f1af5d2fSBarry Smith   PetscTruth          flg;
190dfbe8321SBarry Smith   PetscErrorCode      ierr;
19177431f27SBarry Smith   PetscInt            i;
1922fc52814SBarry Smith   const char          *deft;
193e8105e01SRichard Katz   char                type[256], monfilename[PETSC_MAX_PATH_LEN];
194e8105e01SRichard Katz   PetscViewer         monviewer;
1959b94acceSBarry Smith 
1963a40ed3dSBarry Smith   PetscFunctionBegin;
1974482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
198ca161407SBarry Smith 
199b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
200186905e3SBarry Smith     if (snes->type_name) {
201186905e3SBarry Smith       deft = snes->type_name;
202186905e3SBarry Smith     } else {
2034b27c08aSLois Curfman McInnes       deft = SNESLS;
204d64ed03dSBarry Smith     }
2054bbc92c1SBarry Smith 
206186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
207b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
208d64ed03dSBarry Smith     if (flg) {
209186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
210186905e3SBarry Smith     } else if (!snes->type_name) {
211186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
212d64ed03dSBarry Smith     }
213909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
21493c39befSBarry Smith 
21587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
21670441072SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
217186905e3SBarry Smith 
21887828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
219b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
220b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
22150ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
2225968eb51SBarry Smith     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
2235968eb51SBarry Smith     if (flg) {
2245968eb51SBarry Smith       snes->printreason = PETSC_TRUE;
2255968eb51SBarry Smith     }
226186905e3SBarry Smith 
2277aaed0d8SBarry Smith     ierr = PetscOptionsTruth("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr);
228186905e3SBarry Smith 
229b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
23087828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
23187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
23287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
23387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
23487828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
23587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
236186905e3SBarry Smith 
237b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
23893c39befSBarry Smith     if (flg) {snes->converged = 0;}
239b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
2405cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
241eabae89aSBarry Smith 
242eabae89aSBarry Smith     ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
243e8105e01SRichard Katz     if (flg) {
244e8105e01SRichard Katz       ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr);
245e8105e01SRichard Katz       ierr = SNESSetMonitor(snes,SNESDefaultMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
246e8105e01SRichard Katz     }
247eabae89aSBarry Smith 
248eabae89aSBarry Smith     ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESSetRatioMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
249eabae89aSBarry Smith     if (flg) {
250eabae89aSBarry Smith       ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr);
251eabae89aSBarry Smith       ierr = SNESSetRatioMonitor(snes,monviewer);CHKERRQ(ierr);
252e8105e01SRichard Katz     }
253eabae89aSBarry Smith 
254eabae89aSBarry Smith     ierr = PetscOptionsString("-snes_smonitor","Monitor norm of function (fewer digits)","SNESSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
255eabae89aSBarry Smith     if (flg) {
256eabae89aSBarry Smith       ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr);
257eabae89aSBarry Smith       ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
258eabae89aSBarry Smith     }
259eabae89aSBarry Smith 
260b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
261b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
262b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2637c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
2645ed2d596SBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
2655ed2d596SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
266b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
267186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
268e24b481bSBarry Smith 
269b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2704b27c08aSLois Curfman McInnes     if (flg) {
271186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
272*ae15b995SBarry Smith       ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr);
2739b94acceSBarry Smith     }
274639f9d9dSBarry Smith 
27576b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
27676b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
27776b2cf59SMatthew Knepley     }
27876b2cf59SMatthew Knepley 
279186905e3SBarry Smith     if (snes->setfromoptions) {
280186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
281639f9d9dSBarry Smith     }
282b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2834bbc92c1SBarry Smith 
28494b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
28594b7f48cSBarry Smith   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
28693993e2dSLois Curfman McInnes 
2873a40ed3dSBarry Smith   PetscFunctionReturn(0);
2889b94acceSBarry Smith }
2899b94acceSBarry Smith 
290a847f771SSatish Balay 
2914a2ae208SSatish Balay #undef __FUNCT__
2924a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
2939b94acceSBarry Smith /*@
2949b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2959b94acceSBarry Smith    the nonlinear solvers.
2969b94acceSBarry Smith 
297fee21e36SBarry Smith    Collective on SNES
298fee21e36SBarry Smith 
299c7afd0dbSLois Curfman McInnes    Input Parameters:
300c7afd0dbSLois Curfman McInnes +  snes - the SNES context
301c7afd0dbSLois Curfman McInnes -  usrP - optional user context
302c7afd0dbSLois Curfman McInnes 
30336851e7fSLois Curfman McInnes    Level: intermediate
30436851e7fSLois Curfman McInnes 
3059b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3069b94acceSBarry Smith 
3079b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3089b94acceSBarry Smith @*/
30963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP)
3109b94acceSBarry Smith {
3113a40ed3dSBarry Smith   PetscFunctionBegin;
3124482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3139b94acceSBarry Smith   snes->user		= usrP;
3143a40ed3dSBarry Smith   PetscFunctionReturn(0);
3159b94acceSBarry Smith }
31674679c65SBarry Smith 
3174a2ae208SSatish Balay #undef __FUNCT__
3184a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3199b94acceSBarry Smith /*@C
3209b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3219b94acceSBarry Smith    nonlinear solvers.
3229b94acceSBarry Smith 
323c7afd0dbSLois Curfman McInnes    Not Collective
324c7afd0dbSLois Curfman McInnes 
3259b94acceSBarry Smith    Input Parameter:
3269b94acceSBarry Smith .  snes - SNES context
3279b94acceSBarry Smith 
3289b94acceSBarry Smith    Output Parameter:
3299b94acceSBarry Smith .  usrP - user context
3309b94acceSBarry Smith 
33136851e7fSLois Curfman McInnes    Level: intermediate
33236851e7fSLois Curfman McInnes 
3339b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3349b94acceSBarry Smith 
3359b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3369b94acceSBarry Smith @*/
33763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP)
3389b94acceSBarry Smith {
3393a40ed3dSBarry Smith   PetscFunctionBegin;
3404482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3419b94acceSBarry Smith   *usrP = snes->user;
3423a40ed3dSBarry Smith   PetscFunctionReturn(0);
3439b94acceSBarry Smith }
34474679c65SBarry Smith 
3454a2ae208SSatish Balay #undef __FUNCT__
3464a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3479b94acceSBarry Smith /*@
348c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
349c8228a4eSBarry Smith    at this time.
3509b94acceSBarry Smith 
351c7afd0dbSLois Curfman McInnes    Not Collective
352c7afd0dbSLois Curfman McInnes 
3539b94acceSBarry Smith    Input Parameter:
3549b94acceSBarry Smith .  snes - SNES context
3559b94acceSBarry Smith 
3569b94acceSBarry Smith    Output Parameter:
3579b94acceSBarry Smith .  iter - iteration number
3589b94acceSBarry Smith 
359c8228a4eSBarry Smith    Notes:
360c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
361c8228a4eSBarry Smith 
362c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
36308405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
36408405cd6SLois Curfman McInnes .vb
36508405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
36608405cd6SLois Curfman McInnes       if (!(it % 2)) {
36708405cd6SLois Curfman McInnes         [compute Jacobian here]
36808405cd6SLois Curfman McInnes       }
36908405cd6SLois Curfman McInnes .ve
370c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
37108405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
372c8228a4eSBarry Smith 
37336851e7fSLois Curfman McInnes    Level: intermediate
37436851e7fSLois Curfman McInnes 
3759b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3769b94acceSBarry Smith @*/
37763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter)
3789b94acceSBarry Smith {
3793a40ed3dSBarry Smith   PetscFunctionBegin;
3804482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3814482741eSBarry Smith   PetscValidIntPointer(iter,2);
3829b94acceSBarry Smith   *iter = snes->iter;
3833a40ed3dSBarry Smith   PetscFunctionReturn(0);
3849b94acceSBarry Smith }
38574679c65SBarry Smith 
3864a2ae208SSatish Balay #undef __FUNCT__
3874a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
3889b94acceSBarry Smith /*@
3899b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3909b94acceSBarry Smith    with SNESSSetFunction().
3919b94acceSBarry Smith 
392c7afd0dbSLois Curfman McInnes    Collective on SNES
393c7afd0dbSLois Curfman McInnes 
3949b94acceSBarry Smith    Input Parameter:
3959b94acceSBarry Smith .  snes - SNES context
3969b94acceSBarry Smith 
3979b94acceSBarry Smith    Output Parameter:
3989b94acceSBarry Smith .  fnorm - 2-norm of function
3999b94acceSBarry Smith 
40036851e7fSLois Curfman McInnes    Level: intermediate
40136851e7fSLois Curfman McInnes 
4029b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
403a86d99e1SLois Curfman McInnes 
40408405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
4059b94acceSBarry Smith @*/
40663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
4079b94acceSBarry Smith {
4083a40ed3dSBarry Smith   PetscFunctionBegin;
4094482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4104482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
4119b94acceSBarry Smith   *fnorm = snes->norm;
4123a40ed3dSBarry Smith   PetscFunctionReturn(0);
4139b94acceSBarry Smith }
41474679c65SBarry Smith 
4154a2ae208SSatish Balay #undef __FUNCT__
4164a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
4179b94acceSBarry Smith /*@
4189b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4199b94acceSBarry Smith    attempted by the nonlinear solver.
4209b94acceSBarry Smith 
421c7afd0dbSLois Curfman McInnes    Not Collective
422c7afd0dbSLois Curfman McInnes 
4239b94acceSBarry Smith    Input Parameter:
4249b94acceSBarry Smith .  snes - SNES context
4259b94acceSBarry Smith 
4269b94acceSBarry Smith    Output Parameter:
4279b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4289b94acceSBarry Smith 
429c96a6f78SLois Curfman McInnes    Notes:
430c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
431c96a6f78SLois Curfman McInnes 
43236851e7fSLois Curfman McInnes    Level: intermediate
43336851e7fSLois Curfman McInnes 
4349b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4359b94acceSBarry Smith @*/
43663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* nfails)
4379b94acceSBarry Smith {
4383a40ed3dSBarry Smith   PetscFunctionBegin;
4394482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4404482741eSBarry Smith   PetscValidIntPointer(nfails,2);
44150ffb88aSMatthew Knepley   *nfails = snes->numFailures;
44250ffb88aSMatthew Knepley   PetscFunctionReturn(0);
44350ffb88aSMatthew Knepley }
44450ffb88aSMatthew Knepley 
44550ffb88aSMatthew Knepley #undef __FUNCT__
44650ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
44750ffb88aSMatthew Knepley /*@
44850ffb88aSMatthew Knepley    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
44950ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
45050ffb88aSMatthew Knepley 
45150ffb88aSMatthew Knepley    Not Collective
45250ffb88aSMatthew Knepley 
45350ffb88aSMatthew Knepley    Input Parameters:
45450ffb88aSMatthew Knepley +  snes     - SNES context
45550ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
45650ffb88aSMatthew Knepley 
45750ffb88aSMatthew Knepley    Level: intermediate
45850ffb88aSMatthew Knepley 
45950ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
46050ffb88aSMatthew Knepley @*/
46163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt maxFails)
46250ffb88aSMatthew Knepley {
46350ffb88aSMatthew Knepley   PetscFunctionBegin;
4644482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
46550ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
46650ffb88aSMatthew Knepley   PetscFunctionReturn(0);
46750ffb88aSMatthew Knepley }
46850ffb88aSMatthew Knepley 
46950ffb88aSMatthew Knepley #undef __FUNCT__
47050ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
47150ffb88aSMatthew Knepley /*@
47250ffb88aSMatthew Knepley    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
47350ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
47450ffb88aSMatthew Knepley 
47550ffb88aSMatthew Knepley    Not Collective
47650ffb88aSMatthew Knepley 
47750ffb88aSMatthew Knepley    Input Parameter:
47850ffb88aSMatthew Knepley .  snes     - SNES context
47950ffb88aSMatthew Knepley 
48050ffb88aSMatthew Knepley    Output Parameter:
48150ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
48250ffb88aSMatthew Knepley 
48350ffb88aSMatthew Knepley    Level: intermediate
48450ffb88aSMatthew Knepley 
48550ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
48650ffb88aSMatthew Knepley @*/
48763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *maxFails)
48850ffb88aSMatthew Knepley {
48950ffb88aSMatthew Knepley   PetscFunctionBegin;
4904482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4914482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
49250ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
4933a40ed3dSBarry Smith   PetscFunctionReturn(0);
4949b94acceSBarry Smith }
495a847f771SSatish Balay 
4964a2ae208SSatish Balay #undef __FUNCT__
4974a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations"
498c96a6f78SLois Curfman McInnes /*@
499c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
500c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
501c96a6f78SLois Curfman McInnes 
502c7afd0dbSLois Curfman McInnes    Not Collective
503c7afd0dbSLois Curfman McInnes 
504c96a6f78SLois Curfman McInnes    Input Parameter:
505c96a6f78SLois Curfman McInnes .  snes - SNES context
506c96a6f78SLois Curfman McInnes 
507c96a6f78SLois Curfman McInnes    Output Parameter:
508c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
509c96a6f78SLois Curfman McInnes 
510c96a6f78SLois Curfman McInnes    Notes:
511c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
512c96a6f78SLois Curfman McInnes 
51336851e7fSLois Curfman McInnes    Level: intermediate
51436851e7fSLois Curfman McInnes 
515c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
516c96a6f78SLois Curfman McInnes @*/
51763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberLinearIterations(SNES snes,PetscInt* lits)
518c96a6f78SLois Curfman McInnes {
5193a40ed3dSBarry Smith   PetscFunctionBegin;
5204482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5214482741eSBarry Smith   PetscValidIntPointer(lits,2);
522c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5233a40ed3dSBarry Smith   PetscFunctionReturn(0);
524c96a6f78SLois Curfman McInnes }
525c96a6f78SLois Curfman McInnes 
5264a2ae208SSatish Balay #undef __FUNCT__
52794b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
52852baeb72SSatish Balay /*@
52994b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
5309b94acceSBarry Smith 
53194b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
532c7afd0dbSLois Curfman McInnes 
5339b94acceSBarry Smith    Input Parameter:
5349b94acceSBarry Smith .  snes - the SNES context
5359b94acceSBarry Smith 
5369b94acceSBarry Smith    Output Parameter:
53794b7f48cSBarry Smith .  ksp - the KSP context
5389b94acceSBarry Smith 
5399b94acceSBarry Smith    Notes:
54094b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
5419b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5429b94acceSBarry Smith    KSP and PC contexts as well.
5439b94acceSBarry Smith 
54436851e7fSLois Curfman McInnes    Level: beginner
54536851e7fSLois Curfman McInnes 
54694b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
5479b94acceSBarry Smith 
54894b7f48cSBarry Smith .seealso: KSPGetPC()
5499b94acceSBarry Smith @*/
55063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp)
5519b94acceSBarry Smith {
5523a40ed3dSBarry Smith   PetscFunctionBegin;
5534482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5544482741eSBarry Smith   PetscValidPointer(ksp,2);
55594b7f48cSBarry Smith   *ksp = snes->ksp;
5563a40ed3dSBarry Smith   PetscFunctionReturn(0);
5579b94acceSBarry Smith }
55882bf6240SBarry Smith 
5594a2ae208SSatish Balay #undef __FUNCT__
5604a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
5616849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
562e24b481bSBarry Smith {
563e24b481bSBarry Smith   PetscFunctionBegin;
564e24b481bSBarry Smith   PetscFunctionReturn(0);
565e24b481bSBarry Smith }
566e24b481bSBarry Smith 
5679b94acceSBarry Smith /* -----------------------------------------------------------*/
5684a2ae208SSatish Balay #undef __FUNCT__
5694a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
57052baeb72SSatish Balay /*@
5719b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5729b94acceSBarry Smith 
573c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
574c7afd0dbSLois Curfman McInnes 
575c7afd0dbSLois Curfman McInnes    Input Parameters:
576c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
5779b94acceSBarry Smith 
5789b94acceSBarry Smith    Output Parameter:
5799b94acceSBarry Smith .  outsnes - the new SNES context
5809b94acceSBarry Smith 
581c7afd0dbSLois Curfman McInnes    Options Database Keys:
582c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
583c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
584c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
585c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
586c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
587c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
588c1f60f51SBarry Smith 
58936851e7fSLois Curfman McInnes    Level: beginner
59036851e7fSLois Curfman McInnes 
5919b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5929b94acceSBarry Smith 
5934b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES
5949b94acceSBarry Smith @*/
59563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes)
5969b94acceSBarry Smith {
597dfbe8321SBarry Smith   PetscErrorCode      ierr;
5989b94acceSBarry Smith   SNES                snes;
5999b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
60037fcc0dbSBarry Smith 
6013a40ed3dSBarry Smith   PetscFunctionBegin;
602ed1caa07SMatthew Knepley   PetscValidPointer(outsnes,2);
6038ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
6048ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
6058ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
6068ba1e511SMatthew Knepley #endif
6078ba1e511SMatthew Knepley 
60852e6d16bSBarry Smith   ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
609e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
6109b94acceSBarry Smith   snes->max_its           = 50;
6119750a799SBarry Smith   snes->max_funcs	  = 10000;
6129b94acceSBarry Smith   snes->norm		  = 0.0;
613b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
614b4874afaSBarry Smith   snes->ttol              = 0.0;
61570441072SBarry Smith   snes->abstol		  = 1.e-50;
6169b94acceSBarry Smith   snes->xtol		  = 1.e-8;
6174b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
6189b94acceSBarry Smith   snes->nfuncs            = 0;
61950ffb88aSMatthew Knepley   snes->numFailures       = 0;
62050ffb88aSMatthew Knepley   snes->maxFailures       = 1;
6217a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
622639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6239b94acceSBarry Smith   snes->data              = 0;
6249b94acceSBarry Smith   snes->view              = 0;
6254dc4c822SBarry Smith   snes->setupcalled       = PETSC_FALSE;
626186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
6276f24a144SLois Curfman McInnes   snes->vwork             = 0;
6286f24a144SLois Curfman McInnes   snes->nwork             = 0;
629758f92a0SBarry Smith   snes->conv_hist_len     = 0;
630758f92a0SBarry Smith   snes->conv_hist_max     = 0;
631758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
632758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
633758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
634184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
6359b94acceSBarry Smith 
6369b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
637b0a32e0cSBarry Smith   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
63852e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));CHKERRQ(ierr);
6399b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6409b94acceSBarry Smith   kctx->version     = 2;
6419b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6429b94acceSBarry Smith                              this was too large for some test cases */
6439b94acceSBarry Smith   kctx->rtol_last   = 0;
6449b94acceSBarry Smith   kctx->rtol_max    = .9;
6459b94acceSBarry Smith   kctx->gamma       = 1.0;
6469b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6479b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6489b94acceSBarry Smith   kctx->threshold   = .1;
6499b94acceSBarry Smith   kctx->lresid_last = 0;
6509b94acceSBarry Smith   kctx->norm_last   = 0;
6519b94acceSBarry Smith 
65294b7f48cSBarry Smith   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
65352e6d16bSBarry Smith   ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
6545334005bSBarry Smith 
6559b94acceSBarry Smith   *outsnes = snes;
65600036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
6573a40ed3dSBarry Smith   PetscFunctionReturn(0);
6589b94acceSBarry Smith }
6599b94acceSBarry Smith 
6604a2ae208SSatish Balay #undef __FUNCT__
6614a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
6629b94acceSBarry Smith /*@C
6639b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6649b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6659b94acceSBarry Smith    equations.
6669b94acceSBarry Smith 
667fee21e36SBarry Smith    Collective on SNES
668fee21e36SBarry Smith 
669c7afd0dbSLois Curfman McInnes    Input Parameters:
670c7afd0dbSLois Curfman McInnes +  snes - the SNES context
671c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
672c7afd0dbSLois Curfman McInnes .  r - vector to store function value
673c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
674c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6759b94acceSBarry Smith 
676c7afd0dbSLois Curfman McInnes    Calling sequence of func:
6778d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
678c7afd0dbSLois Curfman McInnes 
679313e4042SLois Curfman McInnes .  f - function vector
680c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6819b94acceSBarry Smith 
6829b94acceSBarry Smith    Notes:
6839b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6849b94acceSBarry Smith $      f'(x) x = -f(x),
685c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6869b94acceSBarry Smith 
68736851e7fSLois Curfman McInnes    Level: beginner
68836851e7fSLois Curfman McInnes 
6899b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6909b94acceSBarry Smith 
691a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6929b94acceSBarry Smith @*/
69363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
6949b94acceSBarry Smith {
6953a40ed3dSBarry Smith   PetscFunctionBegin;
6964482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6974482741eSBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
698c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,r,2);
699184914b5SBarry Smith 
7009b94acceSBarry Smith   snes->computefunction     = func;
7019b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7029b94acceSBarry Smith   snes->funP                = ctx;
7033a40ed3dSBarry Smith   PetscFunctionReturn(0);
7049b94acceSBarry Smith }
7059b94acceSBarry Smith 
7063ab0aad5SBarry Smith /* --------------------------------------------------------------- */
7073ab0aad5SBarry Smith #undef __FUNCT__
7083ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs"
7093ab0aad5SBarry Smith /*@C
7103ab0aad5SBarry Smith    SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set
7113ab0aad5SBarry Smith    it assumes a zero right hand side.
7123ab0aad5SBarry Smith 
7133ab0aad5SBarry Smith    Collective on SNES
7143ab0aad5SBarry Smith 
7153ab0aad5SBarry Smith    Input Parameters:
7163ab0aad5SBarry Smith +  snes - the SNES context
7173ab0aad5SBarry Smith -  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7183ab0aad5SBarry Smith 
7193ab0aad5SBarry Smith    Level: intermediate
7203ab0aad5SBarry Smith 
7213ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side
7223ab0aad5SBarry Smith 
7231096aae1SMatthew Knepley .seealso: SNESGetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7243ab0aad5SBarry Smith @*/
72563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetRhs(SNES snes,Vec rhs)
7263ab0aad5SBarry Smith {
727dfbe8321SBarry Smith   PetscErrorCode ierr;
7283ab0aad5SBarry Smith 
7293ab0aad5SBarry Smith   PetscFunctionBegin;
7303ab0aad5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7313ab0aad5SBarry Smith   if (rhs) {
7323ab0aad5SBarry Smith     PetscValidHeaderSpecific(rhs,VEC_COOKIE,2);
7333ab0aad5SBarry Smith     PetscCheckSameComm(snes,1,rhs,2);
7343ab0aad5SBarry Smith     ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr);
7353ab0aad5SBarry Smith   }
7363ab0aad5SBarry Smith   if (snes->afine) {
7373ab0aad5SBarry Smith     ierr = VecDestroy(snes->afine);CHKERRQ(ierr);
7383ab0aad5SBarry Smith   }
7393ab0aad5SBarry Smith   snes->afine = rhs;
7403ab0aad5SBarry Smith   PetscFunctionReturn(0);
7413ab0aad5SBarry Smith }
7423ab0aad5SBarry Smith 
7434a2ae208SSatish Balay #undef __FUNCT__
7441096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs"
7451096aae1SMatthew Knepley /*@C
7461096aae1SMatthew Knepley    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
7471096aae1SMatthew Knepley    it assumes a zero right hand side.
7481096aae1SMatthew Knepley 
7491096aae1SMatthew Knepley    Collective on SNES
7501096aae1SMatthew Knepley 
7511096aae1SMatthew Knepley    Input Parameter:
7521096aae1SMatthew Knepley .  snes - the SNES context
7531096aae1SMatthew Knepley 
7541096aae1SMatthew Knepley    Output Parameter:
7551096aae1SMatthew Knepley .  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7561096aae1SMatthew Knepley 
7571096aae1SMatthew Knepley    Level: intermediate
7581096aae1SMatthew Knepley 
7591096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side
7601096aae1SMatthew Knepley 
7611096aae1SMatthew Knepley .seealso: SNESSetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7621096aae1SMatthew Knepley @*/
7631096aae1SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs)
7641096aae1SMatthew Knepley {
7651096aae1SMatthew Knepley   PetscFunctionBegin;
7661096aae1SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7671096aae1SMatthew Knepley   PetscValidPointer(rhs,2);
7681096aae1SMatthew Knepley   *rhs = snes->afine;
7691096aae1SMatthew Knepley   PetscFunctionReturn(0);
7701096aae1SMatthew Knepley }
7711096aae1SMatthew Knepley 
7721096aae1SMatthew Knepley #undef __FUNCT__
7734a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
7749b94acceSBarry Smith /*@
77536851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
7769b94acceSBarry Smith                          SNESSetFunction().
7779b94acceSBarry Smith 
778c7afd0dbSLois Curfman McInnes    Collective on SNES
779c7afd0dbSLois Curfman McInnes 
7809b94acceSBarry Smith    Input Parameters:
781c7afd0dbSLois Curfman McInnes +  snes - the SNES context
782c7afd0dbSLois Curfman McInnes -  x - input vector
7839b94acceSBarry Smith 
7849b94acceSBarry Smith    Output Parameter:
7853638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7869b94acceSBarry Smith 
7871bffabb2SLois Curfman McInnes    Notes:
78836851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
78936851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
79036851e7fSLois Curfman McInnes    themselves.
79136851e7fSLois Curfman McInnes 
79236851e7fSLois Curfman McInnes    Level: developer
79336851e7fSLois Curfman McInnes 
7949b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7959b94acceSBarry Smith 
796a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7979b94acceSBarry Smith @*/
79863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y)
7999b94acceSBarry Smith {
800dfbe8321SBarry Smith   PetscErrorCode ierr;
8019b94acceSBarry Smith 
8023a40ed3dSBarry Smith   PetscFunctionBegin;
8034482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8044482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
8054482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
806c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
807c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
808184914b5SBarry Smith 
809d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8101096aae1SMatthew Knepley   if (snes->computefunction) {
811d64ed03dSBarry Smith     PetscStackPush("SNES user function");
812e9a2bbcdSBarry Smith     CHKMEMQ;
81319717074SBarry Smith     ierr = (*snes->computefunction)(snes,x,y,snes->funP);
814e9a2bbcdSBarry Smith     CHKMEMQ;
815d64ed03dSBarry Smith     PetscStackPop;
816d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
81719717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
81819717074SBarry Smith     }
819d5e45103SBarry Smith     CHKERRQ(ierr);
8201096aae1SMatthew Knepley   } else if (snes->afine) {
8211096aae1SMatthew Knepley     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
8221096aae1SMatthew Knepley   } else {
8231096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
8241096aae1SMatthew Knepley   }
8253ab0aad5SBarry Smith   if (snes->afine) {
826016dedfbSBarry Smith     ierr = VecAXPY(y,-1.0,snes->afine);CHKERRQ(ierr);
8273ab0aad5SBarry Smith   }
828ae3c334cSLois Curfman McInnes   snes->nfuncs++;
829d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8303a40ed3dSBarry Smith   PetscFunctionReturn(0);
8319b94acceSBarry Smith }
8329b94acceSBarry Smith 
8334a2ae208SSatish Balay #undef __FUNCT__
8344a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
83562fef451SLois Curfman McInnes /*@
83662fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
83762fef451SLois Curfman McInnes    set with SNESSetJacobian().
83862fef451SLois Curfman McInnes 
839c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
840c7afd0dbSLois Curfman McInnes 
84162fef451SLois Curfman McInnes    Input Parameters:
842c7afd0dbSLois Curfman McInnes +  snes - the SNES context
843c7afd0dbSLois Curfman McInnes -  x - input vector
84462fef451SLois Curfman McInnes 
84562fef451SLois Curfman McInnes    Output Parameters:
846c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
84762fef451SLois Curfman McInnes .  B - optional preconditioning matrix
848c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
849fee21e36SBarry Smith 
85062fef451SLois Curfman McInnes    Notes:
85162fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
85262fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
85362fef451SLois Curfman McInnes 
85494b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
855dc5a77f8SLois Curfman McInnes    flag parameter.
85662fef451SLois Curfman McInnes 
85736851e7fSLois Curfman McInnes    Level: developer
85836851e7fSLois Curfman McInnes 
85962fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
86062fef451SLois Curfman McInnes 
86194b7f48cSBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators()
86262fef451SLois Curfman McInnes @*/
86363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8649b94acceSBarry Smith {
865dfbe8321SBarry Smith   PetscErrorCode ierr;
8663a40ed3dSBarry Smith 
8673a40ed3dSBarry Smith   PetscFunctionBegin;
8684482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8694482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
8704482741eSBarry Smith   PetscValidPointer(flg,5);
871c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
8723a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
873d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
874c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
875d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
876dc67937bSBarry Smith   CHKMEMQ;
8779b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
878dc67937bSBarry Smith   CHKMEMQ;
879d64ed03dSBarry Smith   PetscStackPop;
880d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
8816d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
8824482741eSBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
8834482741eSBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE,4);
8843a40ed3dSBarry Smith   PetscFunctionReturn(0);
8859b94acceSBarry Smith }
8869b94acceSBarry Smith 
8874a2ae208SSatish Balay #undef __FUNCT__
8884a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
8899b94acceSBarry Smith /*@C
8909b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
891044dda88SLois Curfman McInnes    location to store the matrix.
8929b94acceSBarry Smith 
893c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
894c7afd0dbSLois Curfman McInnes 
8959b94acceSBarry Smith    Input Parameters:
896c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8979b94acceSBarry Smith .  A - Jacobian matrix
8989b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8999b94acceSBarry Smith .  func - Jacobian evaluation routine
900c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
9012cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9029b94acceSBarry Smith 
9039b94acceSBarry Smith    Calling sequence of func:
9048d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9059b94acceSBarry Smith 
906c7afd0dbSLois Curfman McInnes +  x - input vector
9079b94acceSBarry Smith .  A - Jacobian matrix
9089b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
909ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
91094b7f48cSBarry Smith    structure (same as flag in KSPSetOperators())
911c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
9129b94acceSBarry Smith 
9139b94acceSBarry Smith    Notes:
91494b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
9152cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
916ac21db08SLois Curfman McInnes 
917ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9189b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9199b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9209b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9219b94acceSBarry Smith    throughout the global iterations.
9229b94acceSBarry Smith 
92336851e7fSLois Curfman McInnes    Level: beginner
92436851e7fSLois Curfman McInnes 
9259b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9269b94acceSBarry Smith 
9273960baedSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor()
9289b94acceSBarry Smith @*/
92963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
9309b94acceSBarry Smith {
931dfbe8321SBarry Smith   PetscErrorCode ierr;
9323a7fca6bSBarry Smith 
9333a40ed3dSBarry Smith   PetscFunctionBegin;
9344482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9354482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
9364482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
937c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
938c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
9393a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
9403a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
9413a7fca6bSBarry Smith   if (A) {
9423a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
9439b94acceSBarry Smith     snes->jacobian = A;
9443a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
9453a7fca6bSBarry Smith   }
9463a7fca6bSBarry Smith   if (B) {
9473a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
9489b94acceSBarry Smith     snes->jacobian_pre = B;
9493a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
9503a7fca6bSBarry Smith   }
95169a612a9SBarry Smith   ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9523a40ed3dSBarry Smith   PetscFunctionReturn(0);
9539b94acceSBarry Smith }
95462fef451SLois Curfman McInnes 
9554a2ae208SSatish Balay #undef __FUNCT__
9564a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
957c2aafc4cSSatish Balay /*@C
958b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
959b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
960b4fd4287SBarry Smith 
961c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
962c7afd0dbSLois Curfman McInnes 
963b4fd4287SBarry Smith    Input Parameter:
964b4fd4287SBarry Smith .  snes - the nonlinear solver context
965b4fd4287SBarry Smith 
966b4fd4287SBarry Smith    Output Parameters:
967c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
968b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
96970e92668SMatthew Knepley .  func - location to put Jacobian function (or PETSC_NULL)
97070e92668SMatthew Knepley -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
971fee21e36SBarry Smith 
97236851e7fSLois Curfman McInnes    Level: advanced
97336851e7fSLois Curfman McInnes 
974b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
975b4fd4287SBarry Smith @*/
97670e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
977b4fd4287SBarry Smith {
9783a40ed3dSBarry Smith   PetscFunctionBegin;
9794482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
980b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
981b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
98200036973SBarry Smith   if (func) *func = snes->computejacobian;
98370e92668SMatthew Knepley   if (ctx)  *ctx  = snes->jacP;
9843a40ed3dSBarry Smith   PetscFunctionReturn(0);
985b4fd4287SBarry Smith }
986b4fd4287SBarry Smith 
9879b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
98863dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
9899b94acceSBarry Smith 
9904a2ae208SSatish Balay #undef __FUNCT__
9914a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
9929b94acceSBarry Smith /*@
9939b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
994272ac6f2SLois Curfman McInnes    of a nonlinear solver.
9959b94acceSBarry Smith 
996fee21e36SBarry Smith    Collective on SNES
997fee21e36SBarry Smith 
998c7afd0dbSLois Curfman McInnes    Input Parameters:
99970e92668SMatthew Knepley .  snes - the SNES context
1000c7afd0dbSLois Curfman McInnes 
1001272ac6f2SLois Curfman McInnes    Notes:
1002272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1003272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1004272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1005272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1006272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1007272ac6f2SLois Curfman McInnes 
100836851e7fSLois Curfman McInnes    Level: advanced
100936851e7fSLois Curfman McInnes 
10109b94acceSBarry Smith .keywords: SNES, nonlinear, setup
10119b94acceSBarry Smith 
10129b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
10139b94acceSBarry Smith @*/
101470e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
10159b94acceSBarry Smith {
1016dfbe8321SBarry Smith   PetscErrorCode ierr;
10174b27c08aSLois Curfman McInnes   PetscTruth     flg, iseqtr;
10183a40ed3dSBarry Smith 
10193a40ed3dSBarry Smith   PetscFunctionBegin;
10204482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10214dc4c822SBarry Smith   if (snes->setupcalled) PetscFunctionReturn(0);
10229b94acceSBarry Smith 
1023b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1024c1f60f51SBarry Smith   /*
1025c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1026dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1027c1f60f51SBarry Smith   */
1028c1f60f51SBarry Smith   if (flg) {
1029c1f60f51SBarry Smith     Mat J;
10305a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10315a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1032*ae15b995SBarry Smith     ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr);
10333a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
10343a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1035c1f60f51SBarry Smith   }
103645fc7adcSBarry Smith 
103703c60df9SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT)
103845fc7adcSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
103945fc7adcSBarry Smith   if (flg) {
104045fc7adcSBarry Smith     Mat J;
104145fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
104245fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
104345fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
104445fc7adcSBarry Smith   }
104532a4b47aSMatthew Knepley #endif
104645fc7adcSBarry Smith 
1047b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1048c1f60f51SBarry Smith   /*
1049dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1050c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1051c1f60f51SBarry Smith    */
105231615d04SLois Curfman McInnes   if (flg) {
1053272ac6f2SLois Curfman McInnes     Mat  J;
1054b5d62d44SBarry Smith     KSP ksp;
105594b7f48cSBarry Smith     PC   pc;
1056f3ef73ceSBarry Smith 
10575a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10583a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1059*ae15b995SBarry Smith     ierr = PetscInfo(snes,"Setting default matrix-free operator and preconditioner routines;\nThat is no preconditioner is being used.\n");CHKERRQ(ierr);
10603a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
10613a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
10623a7fca6bSBarry Smith 
1063f3ef73ceSBarry Smith     /* force no preconditioner */
106494b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1065b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1066a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1067a9815358SBarry Smith     if (!flg) {
1068f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1069272ac6f2SLois Curfman McInnes     }
1070a9815358SBarry Smith   }
1071f3ef73ceSBarry Smith 
10721096aae1SMatthew Knepley   if (!snes->vec_func && !snes->afine) {
10731096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10741096aae1SMatthew Knepley   }
10751096aae1SMatthew Knepley   if (!snes->computefunction && !snes->afine) {
10761096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10771096aae1SMatthew Knepley   }
107829bbc08cSBarry Smith   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1079a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
108029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1081a8c6a408SBarry Smith   }
1082a703fe33SLois Curfman McInnes 
1083a703fe33SLois Curfman McInnes   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
10844b27c08aSLois Curfman McInnes   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
10856831982aSBarry Smith   if (snes->ksp_ewconv && !iseqtr) {
108694b7f48cSBarry Smith     KSP ksp;
108794b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1088186905e3SBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1089a703fe33SLois Curfman McInnes   }
10904b27c08aSLois Curfman McInnes 
1091a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
10927aaed0d8SBarry Smith   snes->setupcalled = PETSC_TRUE;
10933a40ed3dSBarry Smith   PetscFunctionReturn(0);
10949b94acceSBarry Smith }
10959b94acceSBarry Smith 
10964a2ae208SSatish Balay #undef __FUNCT__
10974a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
109852baeb72SSatish Balay /*@
10999b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
11009b94acceSBarry Smith    with SNESCreate().
11019b94acceSBarry Smith 
1102c7afd0dbSLois Curfman McInnes    Collective on SNES
1103c7afd0dbSLois Curfman McInnes 
11049b94acceSBarry Smith    Input Parameter:
11059b94acceSBarry Smith .  snes - the SNES context
11069b94acceSBarry Smith 
110736851e7fSLois Curfman McInnes    Level: beginner
110836851e7fSLois Curfman McInnes 
11099b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11109b94acceSBarry Smith 
111163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11129b94acceSBarry Smith @*/
111363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
11149b94acceSBarry Smith {
11156849ba73SBarry Smith   PetscErrorCode ierr;
11163a40ed3dSBarry Smith 
11173a40ed3dSBarry Smith   PetscFunctionBegin;
11184482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
11193a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1120d4bb536fSBarry Smith 
1121be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
11220f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1123be0abb6dSBarry Smith 
1124e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1125606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
11263a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
11273a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
11283ab0aad5SBarry Smith   if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);}
112994b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1130522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1131d952e501SBarry Smith   ierr = SNESClearMonitor(snes);CHKERRQ(ierr);
1132a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
11333a40ed3dSBarry Smith   PetscFunctionReturn(0);
11349b94acceSBarry Smith }
11359b94acceSBarry Smith 
11369b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11379b94acceSBarry Smith 
11384a2ae208SSatish Balay #undef __FUNCT__
11394a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
11409b94acceSBarry Smith /*@
1141d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11429b94acceSBarry Smith 
1143c7afd0dbSLois Curfman McInnes    Collective on SNES
1144c7afd0dbSLois Curfman McInnes 
11459b94acceSBarry Smith    Input Parameters:
1146c7afd0dbSLois Curfman McInnes +  snes - the SNES context
114770441072SBarry Smith .  abstol - absolute convergence tolerance
114833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
114933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
115033174efeSLois Curfman McInnes            of the change in the solution between steps
115133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1152c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1153fee21e36SBarry Smith 
115433174efeSLois Curfman McInnes    Options Database Keys:
115570441072SBarry Smith +    -snes_atol <abstol> - Sets abstol
1156c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1157c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1158c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1159c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
11609b94acceSBarry Smith 
1161d7a720efSLois Curfman McInnes    Notes:
11629b94acceSBarry Smith    The default maximum number of iterations is 50.
11639b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11649b94acceSBarry Smith 
116536851e7fSLois Curfman McInnes    Level: intermediate
116636851e7fSLois Curfman McInnes 
116733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11689b94acceSBarry Smith 
11692492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance()
11709b94acceSBarry Smith @*/
117163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
11729b94acceSBarry Smith {
11733a40ed3dSBarry Smith   PetscFunctionBegin;
11744482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
117570441072SBarry Smith   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1176d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1177d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1178d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1179d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11803a40ed3dSBarry Smith   PetscFunctionReturn(0);
11819b94acceSBarry Smith }
11829b94acceSBarry Smith 
11834a2ae208SSatish Balay #undef __FUNCT__
11844a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
11859b94acceSBarry Smith /*@
118633174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
118733174efeSLois Curfman McInnes 
1188c7afd0dbSLois Curfman McInnes    Not Collective
1189c7afd0dbSLois Curfman McInnes 
119033174efeSLois Curfman McInnes    Input Parameters:
1191c7afd0dbSLois Curfman McInnes +  snes - the SNES context
119270441072SBarry Smith .  abstol - absolute convergence tolerance
119333174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
119433174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
119533174efeSLois Curfman McInnes            of the change in the solution between steps
119633174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1197c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1198fee21e36SBarry Smith 
119933174efeSLois Curfman McInnes    Notes:
120033174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
120133174efeSLois Curfman McInnes 
120236851e7fSLois Curfman McInnes    Level: intermediate
120336851e7fSLois Curfman McInnes 
120433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
120533174efeSLois Curfman McInnes 
120633174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
120733174efeSLois Curfman McInnes @*/
120863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
120933174efeSLois Curfman McInnes {
12103a40ed3dSBarry Smith   PetscFunctionBegin;
12114482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
121270441072SBarry Smith   if (abstol)  *abstol  = snes->abstol;
121333174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
121433174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
121533174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
121633174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
12173a40ed3dSBarry Smith   PetscFunctionReturn(0);
121833174efeSLois Curfman McInnes }
121933174efeSLois Curfman McInnes 
12204a2ae208SSatish Balay #undef __FUNCT__
12214a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
122233174efeSLois Curfman McInnes /*@
12239b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12249b94acceSBarry Smith 
1225fee21e36SBarry Smith    Collective on SNES
1226fee21e36SBarry Smith 
1227c7afd0dbSLois Curfman McInnes    Input Parameters:
1228c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1229c7afd0dbSLois Curfman McInnes -  tol - tolerance
1230c7afd0dbSLois Curfman McInnes 
12319b94acceSBarry Smith    Options Database Key:
1232c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
12339b94acceSBarry Smith 
123436851e7fSLois Curfman McInnes    Level: intermediate
123536851e7fSLois Curfman McInnes 
12369b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12379b94acceSBarry Smith 
12382492ecdbSBarry Smith .seealso: SNESSetTolerances()
12399b94acceSBarry Smith @*/
124063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
12419b94acceSBarry Smith {
12423a40ed3dSBarry Smith   PetscFunctionBegin;
12434482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12449b94acceSBarry Smith   snes->deltatol = tol;
12453a40ed3dSBarry Smith   PetscFunctionReturn(0);
12469b94acceSBarry Smith }
12479b94acceSBarry Smith 
1248df9fa365SBarry Smith /*
1249df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1250df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1251df9fa365SBarry Smith    macros instead of functions
1252df9fa365SBarry Smith */
12534a2ae208SSatish Balay #undef __FUNCT__
12544a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
125563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1256ce1608b8SBarry Smith {
1257dfbe8321SBarry Smith   PetscErrorCode ierr;
1258ce1608b8SBarry Smith 
1259ce1608b8SBarry Smith   PetscFunctionBegin;
12604482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1261ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1262ce1608b8SBarry Smith   PetscFunctionReturn(0);
1263ce1608b8SBarry Smith }
1264ce1608b8SBarry Smith 
12654a2ae208SSatish Balay #undef __FUNCT__
12664a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
126763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1268df9fa365SBarry Smith {
1269dfbe8321SBarry Smith   PetscErrorCode ierr;
1270df9fa365SBarry Smith 
1271df9fa365SBarry Smith   PetscFunctionBegin;
1272df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1273df9fa365SBarry Smith   PetscFunctionReturn(0);
1274df9fa365SBarry Smith }
1275df9fa365SBarry Smith 
12764a2ae208SSatish Balay #undef __FUNCT__
12774a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
127863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw)
1279df9fa365SBarry Smith {
1280dfbe8321SBarry Smith   PetscErrorCode ierr;
1281df9fa365SBarry Smith 
1282df9fa365SBarry Smith   PetscFunctionBegin;
1283df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1284df9fa365SBarry Smith   PetscFunctionReturn(0);
1285df9fa365SBarry Smith }
1286df9fa365SBarry Smith 
12879b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12889b94acceSBarry Smith 
12894a2ae208SSatish Balay #undef __FUNCT__
12904a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
12919b94acceSBarry Smith /*@C
12925cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
12939b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12949b94acceSBarry Smith    progress.
12959b94acceSBarry Smith 
1296fee21e36SBarry Smith    Collective on SNES
1297fee21e36SBarry Smith 
1298c7afd0dbSLois Curfman McInnes    Input Parameters:
1299c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1300c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1301b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1302e8105e01SRichard Katz           monitor routine (use PETSC_NULL if no context is desired)
1303b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1304b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
13059b94acceSBarry Smith 
1306c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1307a7cc72afSBarry Smith $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1308c7afd0dbSLois Curfman McInnes 
1309c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1310c7afd0dbSLois Curfman McInnes .    its - iteration number
1311c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
131240a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
13139b94acceSBarry Smith 
13149665c990SLois Curfman McInnes    Options Database Keys:
1315c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1316c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1317c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1318c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1319c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1320c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1321c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1322c7afd0dbSLois Curfman McInnes                             the options database.
13239665c990SLois Curfman McInnes 
1324639f9d9dSBarry Smith    Notes:
13256bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13266bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13276bc08f3fSLois Curfman McInnes    order in which they were set.
1328639f9d9dSBarry Smith 
132936851e7fSLois Curfman McInnes    Level: intermediate
133036851e7fSLois Curfman McInnes 
13319b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13329b94acceSBarry Smith 
13335cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
13349b94acceSBarry Smith @*/
133563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
13369b94acceSBarry Smith {
13373a40ed3dSBarry Smith   PetscFunctionBegin;
13384482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1339639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
134029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1341639f9d9dSBarry Smith   }
1342639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1343b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1344639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13453a40ed3dSBarry Smith   PetscFunctionReturn(0);
13469b94acceSBarry Smith }
13479b94acceSBarry Smith 
13484a2ae208SSatish Balay #undef __FUNCT__
13494a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
13505cd90555SBarry Smith /*@C
13515cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
13525cd90555SBarry Smith 
1353c7afd0dbSLois Curfman McInnes    Collective on SNES
1354c7afd0dbSLois Curfman McInnes 
13555cd90555SBarry Smith    Input Parameters:
13565cd90555SBarry Smith .  snes - the SNES context
13575cd90555SBarry Smith 
13581a480d89SAdministrator    Options Database Key:
1359c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1360c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1361c7afd0dbSLois Curfman McInnes     set via the options database
13625cd90555SBarry Smith 
13635cd90555SBarry Smith    Notes:
13645cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
13655cd90555SBarry Smith 
136636851e7fSLois Curfman McInnes    Level: intermediate
136736851e7fSLois Curfman McInnes 
13685cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
13695cd90555SBarry Smith 
13705cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
13715cd90555SBarry Smith @*/
137263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes)
13735cd90555SBarry Smith {
1374d952e501SBarry Smith   PetscErrorCode ierr;
1375d952e501SBarry Smith   PetscInt       i;
1376d952e501SBarry Smith 
13775cd90555SBarry Smith   PetscFunctionBegin;
13784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1379d952e501SBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1380d952e501SBarry Smith     if (snes->monitordestroy[i]) {
1381d952e501SBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1382d952e501SBarry Smith     }
1383d952e501SBarry Smith   }
13845cd90555SBarry Smith   snes->numbermonitors = 0;
13855cd90555SBarry Smith   PetscFunctionReturn(0);
13865cd90555SBarry Smith }
13875cd90555SBarry Smith 
13884a2ae208SSatish Balay #undef __FUNCT__
13894a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
13909b94acceSBarry Smith /*@C
13919b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13929b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13939b94acceSBarry Smith 
1394fee21e36SBarry Smith    Collective on SNES
1395fee21e36SBarry Smith 
1396c7afd0dbSLois Curfman McInnes    Input Parameters:
1397c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1398c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1399c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1400c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
14019b94acceSBarry Smith 
1402c7afd0dbSLois Curfman McInnes    Calling sequence of func:
14036849ba73SBarry Smith $     PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1404c7afd0dbSLois Curfman McInnes 
1405c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1406c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1407184914b5SBarry Smith .    reason - reason for convergence/divergence
1408c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
14094b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
14104b27c08aSLois Curfman McInnes -    f - 2-norm of function
14119b94acceSBarry Smith 
141236851e7fSLois Curfman McInnes    Level: advanced
141336851e7fSLois Curfman McInnes 
14149b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14159b94acceSBarry Smith 
14164b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR()
14179b94acceSBarry Smith @*/
141863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
14199b94acceSBarry Smith {
14203a40ed3dSBarry Smith   PetscFunctionBegin;
14214482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14229b94acceSBarry Smith   (snes)->converged = func;
14239b94acceSBarry Smith   (snes)->cnvP      = cctx;
14243a40ed3dSBarry Smith   PetscFunctionReturn(0);
14259b94acceSBarry Smith }
14269b94acceSBarry Smith 
14274a2ae208SSatish Balay #undef __FUNCT__
14284a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
142952baeb72SSatish Balay /*@
1430184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1431184914b5SBarry Smith 
1432184914b5SBarry Smith    Not Collective
1433184914b5SBarry Smith 
1434184914b5SBarry Smith    Input Parameter:
1435184914b5SBarry Smith .  snes - the SNES context
1436184914b5SBarry Smith 
1437184914b5SBarry Smith    Output Parameter:
1438e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1439184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1440184914b5SBarry Smith 
1441184914b5SBarry Smith    Level: intermediate
1442184914b5SBarry Smith 
1443184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1444184914b5SBarry Smith 
1445184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1446184914b5SBarry Smith 
14474b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1448184914b5SBarry Smith @*/
144963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1450184914b5SBarry Smith {
1451184914b5SBarry Smith   PetscFunctionBegin;
14524482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14534482741eSBarry Smith   PetscValidPointer(reason,2);
1454184914b5SBarry Smith   *reason = snes->reason;
1455184914b5SBarry Smith   PetscFunctionReturn(0);
1456184914b5SBarry Smith }
1457184914b5SBarry Smith 
14584a2ae208SSatish Balay #undef __FUNCT__
14594a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1460c9005455SLois Curfman McInnes /*@
1461c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1462c9005455SLois Curfman McInnes 
1463fee21e36SBarry Smith    Collective on SNES
1464fee21e36SBarry Smith 
1465c7afd0dbSLois Curfman McInnes    Input Parameters:
1466c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1467c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1468cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1469758f92a0SBarry Smith .  na  - size of a and its
147064731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1471758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1472c7afd0dbSLois Curfman McInnes 
1473c9005455SLois Curfman McInnes    Notes:
14744b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1475c9005455SLois Curfman McInnes    at each step.
1476c9005455SLois Curfman McInnes 
1477c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1478c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1479c9005455SLois Curfman McInnes    during the section of code that is being timed.
1480c9005455SLois Curfman McInnes 
148136851e7fSLois Curfman McInnes    Level: intermediate
148236851e7fSLois Curfman McInnes 
1483c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1484758f92a0SBarry Smith 
148508405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1486758f92a0SBarry Smith 
1487c9005455SLois Curfman McInnes @*/
148863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset)
1489c9005455SLois Curfman McInnes {
14903a40ed3dSBarry Smith   PetscFunctionBegin;
14914482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14924482741eSBarry Smith   if (na) PetscValidScalarPointer(a,2);
1493c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1494758f92a0SBarry Smith   snes->conv_hist_its   = its;
1495758f92a0SBarry Smith   snes->conv_hist_max   = na;
1496758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1497758f92a0SBarry Smith   PetscFunctionReturn(0);
1498758f92a0SBarry Smith }
1499758f92a0SBarry Smith 
15004a2ae208SSatish Balay #undef __FUNCT__
15014a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
15020c4c9dddSBarry Smith /*@C
1503758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1504758f92a0SBarry Smith 
1505758f92a0SBarry Smith    Collective on SNES
1506758f92a0SBarry Smith 
1507758f92a0SBarry Smith    Input Parameter:
1508758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1509758f92a0SBarry Smith 
1510758f92a0SBarry Smith    Output Parameters:
1511758f92a0SBarry Smith .  a   - array to hold history
1512758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1513758f92a0SBarry Smith          negative if not converged) for each solve.
1514758f92a0SBarry Smith -  na  - size of a and its
1515758f92a0SBarry Smith 
1516758f92a0SBarry Smith    Notes:
1517758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1518758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1519758f92a0SBarry Smith 
1520758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1521758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1522758f92a0SBarry Smith    during the section of code that is being timed.
1523758f92a0SBarry Smith 
1524758f92a0SBarry Smith    Level: intermediate
1525758f92a0SBarry Smith 
1526758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1527758f92a0SBarry Smith 
1528758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1529758f92a0SBarry Smith 
1530758f92a0SBarry Smith @*/
153163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1532758f92a0SBarry Smith {
1533758f92a0SBarry Smith   PetscFunctionBegin;
15344482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1535758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1536758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1537758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
15383a40ed3dSBarry Smith   PetscFunctionReturn(0);
1539c9005455SLois Curfman McInnes }
1540c9005455SLois Curfman McInnes 
1541e74ef692SMatthew Knepley #undef __FUNCT__
1542e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1543ac226902SBarry Smith /*@C
154476b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
15457e4bb74cSBarry Smith   at the beginning o every iteration of the nonlinear solve. Specifically
15467e4bb74cSBarry Smith   it is called just before the Jacobian is "evaluated".
154776b2cf59SMatthew Knepley 
154876b2cf59SMatthew Knepley   Collective on SNES
154976b2cf59SMatthew Knepley 
155076b2cf59SMatthew Knepley   Input Parameters:
155176b2cf59SMatthew Knepley . snes - The nonlinear solver context
155276b2cf59SMatthew Knepley . func - The function
155376b2cf59SMatthew Knepley 
155476b2cf59SMatthew Knepley   Calling sequence of func:
1555b5d30489SBarry Smith . func (SNES snes, PetscInt step);
155676b2cf59SMatthew Knepley 
155776b2cf59SMatthew Knepley . step - The current step of the iteration
155876b2cf59SMatthew Knepley 
155976b2cf59SMatthew Knepley   Level: intermediate
156076b2cf59SMatthew Knepley 
156176b2cf59SMatthew Knepley .keywords: SNES, update
1562b5d30489SBarry Smith 
156376b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
156476b2cf59SMatthew Knepley @*/
156563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
156676b2cf59SMatthew Knepley {
156776b2cf59SMatthew Knepley   PetscFunctionBegin;
15684482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
156976b2cf59SMatthew Knepley   snes->update = func;
157076b2cf59SMatthew Knepley   PetscFunctionReturn(0);
157176b2cf59SMatthew Knepley }
157276b2cf59SMatthew Knepley 
1573e74ef692SMatthew Knepley #undef __FUNCT__
1574e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
157576b2cf59SMatthew Knepley /*@
157676b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
157776b2cf59SMatthew Knepley 
157876b2cf59SMatthew Knepley   Not collective
157976b2cf59SMatthew Knepley 
158076b2cf59SMatthew Knepley   Input Parameters:
158176b2cf59SMatthew Knepley . snes - The nonlinear solver context
158276b2cf59SMatthew Knepley . step - The current step of the iteration
158376b2cf59SMatthew Knepley 
1584205452f4SMatthew Knepley   Level: intermediate
1585205452f4SMatthew Knepley 
158676b2cf59SMatthew Knepley .keywords: SNES, update
158776b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
158876b2cf59SMatthew Knepley @*/
158963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
159076b2cf59SMatthew Knepley {
159176b2cf59SMatthew Knepley   PetscFunctionBegin;
159276b2cf59SMatthew Knepley   PetscFunctionReturn(0);
159376b2cf59SMatthew Knepley }
159476b2cf59SMatthew Knepley 
15954a2ae208SSatish Balay #undef __FUNCT__
15964a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
15979b94acceSBarry Smith /*
15989b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15999b94acceSBarry Smith    positive parameter delta.
16009b94acceSBarry Smith 
16019b94acceSBarry Smith     Input Parameters:
1602c7afd0dbSLois Curfman McInnes +   snes - the SNES context
16039b94acceSBarry Smith .   y - approximate solution of linear system
16049b94acceSBarry Smith .   fnorm - 2-norm of current function
1605c7afd0dbSLois Curfman McInnes -   delta - trust region size
16069b94acceSBarry Smith 
16079b94acceSBarry Smith     Output Parameters:
1608c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
16099b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
16109b94acceSBarry Smith     region, and exceeds zero otherwise.
1611c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
16129b94acceSBarry Smith 
16139b94acceSBarry Smith     Note:
16144b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
16159b94acceSBarry Smith     is set to be the maximum allowable step size.
16169b94acceSBarry Smith 
16179b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
16189b94acceSBarry Smith */
1619dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
16209b94acceSBarry Smith {
1621064f8208SBarry Smith   PetscReal      nrm;
1622ea709b57SSatish Balay   PetscScalar    cnorm;
1623dfbe8321SBarry Smith   PetscErrorCode ierr;
16243a40ed3dSBarry Smith 
16253a40ed3dSBarry Smith   PetscFunctionBegin;
16264482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16274482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1628c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1629184914b5SBarry Smith 
1630064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1631064f8208SBarry Smith   if (nrm > *delta) {
1632064f8208SBarry Smith      nrm = *delta/nrm;
1633064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1634064f8208SBarry Smith      cnorm = nrm;
16352dcb1b2aSMatthew Knepley      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
16369b94acceSBarry Smith      *ynorm = *delta;
16379b94acceSBarry Smith   } else {
16389b94acceSBarry Smith      *gpnorm = 0.0;
1639064f8208SBarry Smith      *ynorm = nrm;
16409b94acceSBarry Smith   }
16413a40ed3dSBarry Smith   PetscFunctionReturn(0);
16429b94acceSBarry Smith }
16439b94acceSBarry Smith 
16444a2ae208SSatish Balay #undef __FUNCT__
16454a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
16469b94acceSBarry Smith /*@
1647f69a0ea3SMatthew Knepley    SNESSolve - Solves a nonlinear system F(x) = b.
1648f69a0ea3SMatthew Knepley    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
16499b94acceSBarry Smith 
1650c7afd0dbSLois Curfman McInnes    Collective on SNES
1651c7afd0dbSLois Curfman McInnes 
1652b2002411SLois Curfman McInnes    Input Parameters:
1653c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1654f69a0ea3SMatthew Knepley .  b - the constant part of the equation, or PETSC_NULL to use zero.
165570e92668SMatthew Knepley -  x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution()
16569b94acceSBarry Smith 
1657b2002411SLois Curfman McInnes    Notes:
16588ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
16598ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16608ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16618ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16628ddd3da0SLois Curfman McInnes 
166336851e7fSLois Curfman McInnes    Level: beginner
166436851e7fSLois Curfman McInnes 
16659b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16669b94acceSBarry Smith 
166770e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution()
16689b94acceSBarry Smith @*/
1669f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
16709b94acceSBarry Smith {
1671dfbe8321SBarry Smith   PetscErrorCode ierr;
1672f1af5d2fSBarry Smith   PetscTruth     flg;
1673eabae89aSBarry Smith   char           filename[PETSC_MAX_PATH_LEN];
1674eabae89aSBarry Smith   PetscViewer    viewer;
1675052efed2SBarry Smith 
16763a40ed3dSBarry Smith   PetscFunctionBegin;
16774482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16781302d50aSBarry Smith   if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
167974637425SBarry Smith 
1680f69a0ea3SMatthew Knepley   if (b) {
1681f69a0ea3SMatthew Knepley     ierr = SNESSetRhs(snes, b); CHKERRQ(ierr);
16821096aae1SMatthew Knepley     if (!snes->vec_func) {
16831096aae1SMatthew Knepley       Vec r;
16841096aae1SMatthew Knepley 
16851096aae1SMatthew Knepley       ierr = VecDuplicate(b, &r); CHKERRQ(ierr);
16861096aae1SMatthew Knepley       ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr);
16871096aae1SMatthew Knepley     }
1688f69a0ea3SMatthew Knepley   }
168970e92668SMatthew Knepley   if (x) {
1690f69a0ea3SMatthew Knepley     PetscValidHeaderSpecific(x,VEC_COOKIE,3);
1691f69a0ea3SMatthew Knepley     PetscCheckSameComm(snes,1,x,3);
169270e92668SMatthew Knepley   } else {
169370e92668SMatthew Knepley     ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr);
169470e92668SMatthew Knepley     if (!x) {
169570e92668SMatthew Knepley       ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr);
169670e92668SMatthew Knepley     }
169770e92668SMatthew Knepley   }
169870e92668SMatthew Knepley   snes->vec_sol = snes->vec_sol_always = x;
169970e92668SMatthew Knepley   if (!snes->setupcalled) {
170070e92668SMatthew Knepley     ierr = SNESSetUp(snes);CHKERRQ(ierr);
170170e92668SMatthew Knepley   }
1702abc0a331SBarry Smith   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
1703d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
170450ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1705d5e45103SBarry Smith 
1706d5e45103SBarry Smith   ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN);
1707d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
170838f152feSBarry Smith     /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */
170919717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr);
1710d5e45103SBarry Smith   } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) {
1711d5e45103SBarry Smith     /* translate exception into SNES not converged reason */
1712d5e45103SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
171338f152feSBarry Smith     ierr = 0;
171438f152feSBarry Smith   }
171538f152feSBarry Smith   CHKERRQ(ierr);
1716d5e45103SBarry Smith 
1717d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1718eabae89aSBarry Smith   ierr = PetscOptionsGetString(snes->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1719eabae89aSBarry Smith   if (flg && !PetscPreLoadingOn) {
1720eabae89aSBarry Smith     ierr = PetscViewerASCIIOpen(snes->comm,filename,&viewer);CHKERRQ(ierr);
1721eabae89aSBarry Smith     ierr = SNESView(snes,viewer);CHKERRQ(ierr);
1722eabae89aSBarry Smith     ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr);
1723eabae89aSBarry Smith   }
1724eabae89aSBarry Smith 
1725da9b6338SBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1726da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
17275968eb51SBarry Smith   if (snes->printreason) {
17285968eb51SBarry Smith     if (snes->reason > 0) {
17299dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17305968eb51SBarry Smith     } else {
17319dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17325968eb51SBarry Smith     }
17335968eb51SBarry Smith   }
17345968eb51SBarry Smith 
17353a40ed3dSBarry Smith   PetscFunctionReturn(0);
17369b94acceSBarry Smith }
17379b94acceSBarry Smith 
17389b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
17399b94acceSBarry Smith 
17404a2ae208SSatish Balay #undef __FUNCT__
17414a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
174282bf6240SBarry Smith /*@C
17434b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
17449b94acceSBarry Smith 
1745fee21e36SBarry Smith    Collective on SNES
1746fee21e36SBarry Smith 
1747c7afd0dbSLois Curfman McInnes    Input Parameters:
1748c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1749454a90a3SBarry Smith -  type - a known method
1750c7afd0dbSLois Curfman McInnes 
1751c7afd0dbSLois Curfman McInnes    Options Database Key:
1752454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1753c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1754ae12b187SLois Curfman McInnes 
17559b94acceSBarry Smith    Notes:
1756e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
17574b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
1758c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17594b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
1760c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17619b94acceSBarry Smith 
1762ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1763ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1764ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1765ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1766ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1767ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1768ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1769ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1770ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1771b0a32e0cSBarry Smith   appropriate method.
177236851e7fSLois Curfman McInnes 
177336851e7fSLois Curfman McInnes   Level: intermediate
1774a703fe33SLois Curfman McInnes 
1775454a90a3SBarry Smith .keywords: SNES, set, type
1776435da068SBarry Smith 
1777435da068SBarry Smith .seealso: SNESType, SNESCreate()
1778435da068SBarry Smith 
17799b94acceSBarry Smith @*/
1780e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type)
17819b94acceSBarry Smith {
1782dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
17836831982aSBarry Smith   PetscTruth     match;
17843a40ed3dSBarry Smith 
17853a40ed3dSBarry Smith   PetscFunctionBegin;
17864482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17874482741eSBarry Smith   PetscValidCharPointer(type,2);
178882bf6240SBarry Smith 
17896831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
17900f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
179182bf6240SBarry Smith 
179282bf6240SBarry Smith   if (snes->setupcalled) {
17934dc4c822SBarry Smith     snes->setupcalled = PETSC_FALSE;
1794e1311b90SBarry Smith     ierr              = (*(snes)->destroy)(snes);CHKERRQ(ierr);
179582bf6240SBarry Smith     snes->data        = 0;
179682bf6240SBarry Smith   }
17977d1a2b2bSBarry Smith 
17989b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
179982bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1800b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1801958c9bccSBarry Smith   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
1802606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
180382bf6240SBarry Smith   snes->data = 0;
18043a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1805454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
18063a40ed3dSBarry Smith   PetscFunctionReturn(0);
18079b94acceSBarry Smith }
18089b94acceSBarry Smith 
1809a847f771SSatish Balay 
18109b94acceSBarry Smith /* --------------------------------------------------------------------- */
18114a2ae208SSatish Balay #undef __FUNCT__
18124a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
181352baeb72SSatish Balay /*@
18149b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1815f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
18169b94acceSBarry Smith 
1817fee21e36SBarry Smith    Not Collective
1818fee21e36SBarry Smith 
181936851e7fSLois Curfman McInnes    Level: advanced
182036851e7fSLois Curfman McInnes 
18219b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
18229b94acceSBarry Smith 
18239b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
18249b94acceSBarry Smith @*/
182563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
18269b94acceSBarry Smith {
1827dfbe8321SBarry Smith   PetscErrorCode ierr;
182882bf6240SBarry Smith 
18293a40ed3dSBarry Smith   PetscFunctionBegin;
183082bf6240SBarry Smith   if (SNESList) {
1831b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
183282bf6240SBarry Smith     SNESList = 0;
18339b94acceSBarry Smith   }
18344c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
18353a40ed3dSBarry Smith   PetscFunctionReturn(0);
18369b94acceSBarry Smith }
18379b94acceSBarry Smith 
18384a2ae208SSatish Balay #undef __FUNCT__
18394a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
18409b94acceSBarry Smith /*@C
18419a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
18429b94acceSBarry Smith 
1843c7afd0dbSLois Curfman McInnes    Not Collective
1844c7afd0dbSLois Curfman McInnes 
18459b94acceSBarry Smith    Input Parameter:
18464b0e389bSBarry Smith .  snes - nonlinear solver context
18479b94acceSBarry Smith 
18489b94acceSBarry Smith    Output Parameter:
18493a7fca6bSBarry Smith .  type - SNES method (a character string)
18509b94acceSBarry Smith 
185136851e7fSLois Curfman McInnes    Level: intermediate
185236851e7fSLois Curfman McInnes 
1853454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
18549b94acceSBarry Smith @*/
185563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type)
18569b94acceSBarry Smith {
18573a40ed3dSBarry Smith   PetscFunctionBegin;
18584482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18594482741eSBarry Smith   PetscValidPointer(type,2);
1860454a90a3SBarry Smith   *type = snes->type_name;
18613a40ed3dSBarry Smith   PetscFunctionReturn(0);
18629b94acceSBarry Smith }
18639b94acceSBarry Smith 
18644a2ae208SSatish Balay #undef __FUNCT__
18654a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
186652baeb72SSatish Balay /*@
18679b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
18689b94acceSBarry Smith    stored.
18699b94acceSBarry Smith 
1870c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1871c7afd0dbSLois Curfman McInnes 
18729b94acceSBarry Smith    Input Parameter:
18739b94acceSBarry Smith .  snes - the SNES context
18749b94acceSBarry Smith 
18759b94acceSBarry Smith    Output Parameter:
18769b94acceSBarry Smith .  x - the solution
18779b94acceSBarry Smith 
187870e92668SMatthew Knepley    Level: intermediate
187936851e7fSLois Curfman McInnes 
18809b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18819b94acceSBarry Smith 
188270e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
18839b94acceSBarry Smith @*/
188463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
18859b94acceSBarry Smith {
18863a40ed3dSBarry Smith   PetscFunctionBegin;
18874482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18884482741eSBarry Smith   PetscValidPointer(x,2);
18899b94acceSBarry Smith   *x = snes->vec_sol_always;
18903a40ed3dSBarry Smith   PetscFunctionReturn(0);
18919b94acceSBarry Smith }
18929b94acceSBarry Smith 
18934a2ae208SSatish Balay #undef __FUNCT__
189470e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution"
18957e4bb74cSBarry Smith /*@
189670e92668SMatthew Knepley    SNESSetSolution - Sets the vector where the approximate solution is stored.
189770e92668SMatthew Knepley 
189870e92668SMatthew Knepley    Not Collective, but Vec is parallel if SNES is parallel
189970e92668SMatthew Knepley 
190070e92668SMatthew Knepley    Input Parameters:
190170e92668SMatthew Knepley +  snes - the SNES context
190270e92668SMatthew Knepley -  x - the solution
190370e92668SMatthew Knepley 
190470e92668SMatthew Knepley    Output Parameter:
190570e92668SMatthew Knepley 
190670e92668SMatthew Knepley    Level: intermediate
190770e92668SMatthew Knepley 
190842219521SBarry Smith    Notes: this is not normally used, rather one simply calls SNESSolve() with
190942219521SBarry Smith           the appropriate solution vector.
191042219521SBarry Smith 
191170e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution
191270e92668SMatthew Knepley 
191370e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
191470e92668SMatthew Knepley @*/
191570e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x)
191670e92668SMatthew Knepley {
191770e92668SMatthew Knepley   PetscFunctionBegin;
191870e92668SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
191970e92668SMatthew Knepley   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
192070e92668SMatthew Knepley   PetscCheckSameComm(snes,1,x,2);
192170e92668SMatthew Knepley   snes->vec_sol_always = x;
192270e92668SMatthew Knepley   PetscFunctionReturn(0);
192370e92668SMatthew Knepley }
192470e92668SMatthew Knepley 
192570e92668SMatthew Knepley #undef __FUNCT__
19264a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
192752baeb72SSatish Balay /*@
19289b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
19299b94acceSBarry Smith    stored.
19309b94acceSBarry Smith 
1931c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1932c7afd0dbSLois Curfman McInnes 
19339b94acceSBarry Smith    Input Parameter:
19349b94acceSBarry Smith .  snes - the SNES context
19359b94acceSBarry Smith 
19369b94acceSBarry Smith    Output Parameter:
19379b94acceSBarry Smith .  x - the solution update
19389b94acceSBarry Smith 
193936851e7fSLois Curfman McInnes    Level: advanced
194036851e7fSLois Curfman McInnes 
19419b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
19429b94acceSBarry Smith 
19439b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
19449b94acceSBarry Smith @*/
194563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
19469b94acceSBarry Smith {
19473a40ed3dSBarry Smith   PetscFunctionBegin;
19484482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19494482741eSBarry Smith   PetscValidPointer(x,2);
19509b94acceSBarry Smith   *x = snes->vec_sol_update_always;
19513a40ed3dSBarry Smith   PetscFunctionReturn(0);
19529b94acceSBarry Smith }
19539b94acceSBarry Smith 
19544a2ae208SSatish Balay #undef __FUNCT__
19554a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
19569b94acceSBarry Smith /*@C
19573638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19589b94acceSBarry Smith 
1959c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1960c7afd0dbSLois Curfman McInnes 
19619b94acceSBarry Smith    Input Parameter:
19629b94acceSBarry Smith .  snes - the SNES context
19639b94acceSBarry Smith 
19649b94acceSBarry Smith    Output Parameter:
19657bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
196670e92668SMatthew Knepley .  func - the function (or PETSC_NULL)
196770e92668SMatthew Knepley -  ctx - the function context (or PETSC_NULL)
19689b94acceSBarry Smith 
196936851e7fSLois Curfman McInnes    Level: advanced
197036851e7fSLois Curfman McInnes 
1971a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19729b94acceSBarry Smith 
19734b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
19749b94acceSBarry Smith @*/
197570e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
19769b94acceSBarry Smith {
19773a40ed3dSBarry Smith   PetscFunctionBegin;
19784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19797bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
198000036973SBarry Smith   if (func) *func = snes->computefunction;
198170e92668SMatthew Knepley   if (ctx)  *ctx  = snes->funP;
19823a40ed3dSBarry Smith   PetscFunctionReturn(0);
19839b94acceSBarry Smith }
19849b94acceSBarry Smith 
19854a2ae208SSatish Balay #undef __FUNCT__
19864a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
19873c7409f5SSatish Balay /*@C
19883c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1989d850072dSLois Curfman McInnes    SNES options in the database.
19903c7409f5SSatish Balay 
1991fee21e36SBarry Smith    Collective on SNES
1992fee21e36SBarry Smith 
1993c7afd0dbSLois Curfman McInnes    Input Parameter:
1994c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1995c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
1996c7afd0dbSLois Curfman McInnes 
1997d850072dSLois Curfman McInnes    Notes:
1998a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1999c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2000d850072dSLois Curfman McInnes 
200136851e7fSLois Curfman McInnes    Level: advanced
200236851e7fSLois Curfman McInnes 
20033c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2004a86d99e1SLois Curfman McInnes 
2005a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
20063c7409f5SSatish Balay @*/
200763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
20083c7409f5SSatish Balay {
2009dfbe8321SBarry Smith   PetscErrorCode ierr;
20103c7409f5SSatish Balay 
20113a40ed3dSBarry Smith   PetscFunctionBegin;
20124482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2013639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
201494b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20153a40ed3dSBarry Smith   PetscFunctionReturn(0);
20163c7409f5SSatish Balay }
20173c7409f5SSatish Balay 
20184a2ae208SSatish Balay #undef __FUNCT__
20194a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
20203c7409f5SSatish Balay /*@C
2021f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2022d850072dSLois Curfman McInnes    SNES options in the database.
20233c7409f5SSatish Balay 
2024fee21e36SBarry Smith    Collective on SNES
2025fee21e36SBarry Smith 
2026c7afd0dbSLois Curfman McInnes    Input Parameters:
2027c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2028c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2029c7afd0dbSLois Curfman McInnes 
2030d850072dSLois Curfman McInnes    Notes:
2031a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2032c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2033d850072dSLois Curfman McInnes 
203436851e7fSLois Curfman McInnes    Level: advanced
203536851e7fSLois Curfman McInnes 
20363c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2037a86d99e1SLois Curfman McInnes 
2038a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20393c7409f5SSatish Balay @*/
204063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
20413c7409f5SSatish Balay {
2042dfbe8321SBarry Smith   PetscErrorCode ierr;
20433c7409f5SSatish Balay 
20443a40ed3dSBarry Smith   PetscFunctionBegin;
20454482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2046639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
204794b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20483a40ed3dSBarry Smith   PetscFunctionReturn(0);
20493c7409f5SSatish Balay }
20503c7409f5SSatish Balay 
20514a2ae208SSatish Balay #undef __FUNCT__
20524a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
20539ab63eb5SSatish Balay /*@C
20543c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20553c7409f5SSatish Balay    SNES options in the database.
20563c7409f5SSatish Balay 
2057c7afd0dbSLois Curfman McInnes    Not Collective
2058c7afd0dbSLois Curfman McInnes 
20593c7409f5SSatish Balay    Input Parameter:
20603c7409f5SSatish Balay .  snes - the SNES context
20613c7409f5SSatish Balay 
20623c7409f5SSatish Balay    Output Parameter:
20633c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20643c7409f5SSatish Balay 
20659ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
20669ab63eb5SSatish Balay    sufficient length to hold the prefix.
20679ab63eb5SSatish Balay 
206836851e7fSLois Curfman McInnes    Level: advanced
206936851e7fSLois Curfman McInnes 
20703c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2071a86d99e1SLois Curfman McInnes 
2072a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20733c7409f5SSatish Balay @*/
2074e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
20753c7409f5SSatish Balay {
2076dfbe8321SBarry Smith   PetscErrorCode ierr;
20773c7409f5SSatish Balay 
20783a40ed3dSBarry Smith   PetscFunctionBegin;
20794482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2080639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
20813a40ed3dSBarry Smith   PetscFunctionReturn(0);
20823c7409f5SSatish Balay }
20833c7409f5SSatish Balay 
2084b2002411SLois Curfman McInnes 
20854a2ae208SSatish Balay #undef __FUNCT__
20864a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
20873cea93caSBarry Smith /*@C
20883cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
20893cea93caSBarry Smith 
20907f6c08e0SMatthew Knepley   Level: advanced
20913cea93caSBarry Smith @*/
209263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2093b2002411SLois Curfman McInnes {
2094e2d1d2b7SBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
2095dfbe8321SBarry Smith   PetscErrorCode ierr;
2096b2002411SLois Curfman McInnes 
2097b2002411SLois Curfman McInnes   PetscFunctionBegin;
2098b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2099c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2100b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2101b2002411SLois Curfman McInnes }
2102da9b6338SBarry Smith 
2103da9b6338SBarry Smith #undef __FUNCT__
2104da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
210563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2106da9b6338SBarry Smith {
2107dfbe8321SBarry Smith   PetscErrorCode ierr;
210877431f27SBarry Smith   PetscInt       N,i,j;
2109da9b6338SBarry Smith   Vec            u,uh,fh;
2110da9b6338SBarry Smith   PetscScalar    value;
2111da9b6338SBarry Smith   PetscReal      norm;
2112da9b6338SBarry Smith 
2113da9b6338SBarry Smith   PetscFunctionBegin;
2114da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2115da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2116da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2117da9b6338SBarry Smith 
2118da9b6338SBarry Smith   /* currently only works for sequential */
2119da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2120da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2121da9b6338SBarry Smith   for (i=0; i<N; i++) {
2122da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
212377431f27SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2124da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2125ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2126da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
21273ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2128da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
212977431f27SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2130da9b6338SBarry Smith       value = -value;
2131da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2132da9b6338SBarry Smith     }
2133da9b6338SBarry Smith   }
2134da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2135da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2136da9b6338SBarry Smith   PetscFunctionReturn(0);
2137da9b6338SBarry Smith }
2138