xref: /petsc/src/snes/interface/snes.c (revision e8105e01a41fc5d9aa83727d0e2b1a76700e3063)
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);
7777d8c4bbSBarry 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);
85b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
86b0a32e0cSBarry 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
156*e8105e01SRichard Katz .  -snes_monitor <optional filename> - prints residual norm at each iteration. if no
157*e8105e01SRichard 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;
193*e8105e01SRichard Katz   char                type[256], monfilename[PETSC_MAX_PATH_LEN];
194*e8105e01SRichard Katz   PetscViewer         monviewer;
195*e8105e01SRichard Katz   size_t              len;
1969b94acceSBarry Smith 
1973a40ed3dSBarry Smith   PetscFunctionBegin;
1984482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
199ca161407SBarry Smith 
200b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
201186905e3SBarry Smith     if (snes->type_name) {
202186905e3SBarry Smith       deft = snes->type_name;
203186905e3SBarry Smith     } else {
2044b27c08aSLois Curfman McInnes       deft = SNESLS;
205d64ed03dSBarry Smith     }
2064bbc92c1SBarry Smith 
207186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
208b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
209d64ed03dSBarry Smith     if (flg) {
210186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
211186905e3SBarry Smith     } else if (!snes->type_name) {
212186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
213d64ed03dSBarry Smith     }
214909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
21593c39befSBarry Smith 
21687828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
21770441072SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
218186905e3SBarry Smith 
21987828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
220b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
221b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
22250ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
2235968eb51SBarry Smith     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
2245968eb51SBarry Smith     if (flg) {
2255968eb51SBarry Smith       snes->printreason = PETSC_TRUE;
2265968eb51SBarry Smith     }
227186905e3SBarry Smith 
2287aaed0d8SBarry 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);
229186905e3SBarry Smith 
230b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
23187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
23287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
23387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
23487828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
23587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
23687828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
237186905e3SBarry Smith 
238b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
23993c39befSBarry Smith     if (flg) {snes->converged = 0;}
240b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
2415cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
242b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
243*e8105e01SRichard Katz     if (flg) {
244*e8105e01SRichard Katz       ierr = PetscOptionsGetString(PETSC_NULL,"-snes_monitor",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
245*e8105e01SRichard Katz       ierr = PetscStrlen(monfilename,&len);CHKERRQ(ierr);
246*e8105e01SRichard Katz       if (len>0) {
247*e8105e01SRichard Katz 	ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr);
248*e8105e01SRichard Katz 	ierr = SNESSetMonitor(snes,SNESDefaultMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr);
249*e8105e01SRichard Katz       } else {
250*e8105e01SRichard Katz 	ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);
251*e8105e01SRichard Katz       }
252*e8105e01SRichard Katz     }
2533a7fca6bSBarry Smith     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
2543a7fca6bSBarry Smith     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
255b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
256b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
257b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
258b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
259b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2607c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
2615ed2d596SBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
2625ed2d596SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
263b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
264186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
265e24b481bSBarry Smith 
266b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2674b27c08aSLois Curfman McInnes     if (flg) {
268186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
26963ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"));CHKERRQ(ierr);
2709b94acceSBarry Smith     }
271639f9d9dSBarry Smith 
27276b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
27376b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
27476b2cf59SMatthew Knepley     }
27576b2cf59SMatthew Knepley 
276186905e3SBarry Smith     if (snes->setfromoptions) {
277186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
278639f9d9dSBarry Smith     }
279639f9d9dSBarry Smith 
280b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2814bbc92c1SBarry Smith 
28294b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
28394b7f48cSBarry Smith   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
28493993e2dSLois Curfman McInnes 
2853a40ed3dSBarry Smith   PetscFunctionReturn(0);
2869b94acceSBarry Smith }
2879b94acceSBarry Smith 
288a847f771SSatish Balay 
2894a2ae208SSatish Balay #undef __FUNCT__
2904a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
2919b94acceSBarry Smith /*@
2929b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2939b94acceSBarry Smith    the nonlinear solvers.
2949b94acceSBarry Smith 
295fee21e36SBarry Smith    Collective on SNES
296fee21e36SBarry Smith 
297c7afd0dbSLois Curfman McInnes    Input Parameters:
298c7afd0dbSLois Curfman McInnes +  snes - the SNES context
299c7afd0dbSLois Curfman McInnes -  usrP - optional user context
300c7afd0dbSLois Curfman McInnes 
30136851e7fSLois Curfman McInnes    Level: intermediate
30236851e7fSLois Curfman McInnes 
3039b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3049b94acceSBarry Smith 
3059b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3069b94acceSBarry Smith @*/
30763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP)
3089b94acceSBarry Smith {
3093a40ed3dSBarry Smith   PetscFunctionBegin;
3104482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3119b94acceSBarry Smith   snes->user		= usrP;
3123a40ed3dSBarry Smith   PetscFunctionReturn(0);
3139b94acceSBarry Smith }
31474679c65SBarry Smith 
3154a2ae208SSatish Balay #undef __FUNCT__
3164a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3179b94acceSBarry Smith /*@C
3189b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3199b94acceSBarry Smith    nonlinear solvers.
3209b94acceSBarry Smith 
321c7afd0dbSLois Curfman McInnes    Not Collective
322c7afd0dbSLois Curfman McInnes 
3239b94acceSBarry Smith    Input Parameter:
3249b94acceSBarry Smith .  snes - SNES context
3259b94acceSBarry Smith 
3269b94acceSBarry Smith    Output Parameter:
3279b94acceSBarry Smith .  usrP - user context
3289b94acceSBarry Smith 
32936851e7fSLois Curfman McInnes    Level: intermediate
33036851e7fSLois Curfman McInnes 
3319b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3329b94acceSBarry Smith 
3339b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3349b94acceSBarry Smith @*/
33563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP)
3369b94acceSBarry Smith {
3373a40ed3dSBarry Smith   PetscFunctionBegin;
3384482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3399b94acceSBarry Smith   *usrP = snes->user;
3403a40ed3dSBarry Smith   PetscFunctionReturn(0);
3419b94acceSBarry Smith }
34274679c65SBarry Smith 
3434a2ae208SSatish Balay #undef __FUNCT__
3444a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3459b94acceSBarry Smith /*@
346c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
347c8228a4eSBarry Smith    at this time.
3489b94acceSBarry Smith 
349c7afd0dbSLois Curfman McInnes    Not Collective
350c7afd0dbSLois Curfman McInnes 
3519b94acceSBarry Smith    Input Parameter:
3529b94acceSBarry Smith .  snes - SNES context
3539b94acceSBarry Smith 
3549b94acceSBarry Smith    Output Parameter:
3559b94acceSBarry Smith .  iter - iteration number
3569b94acceSBarry Smith 
357c8228a4eSBarry Smith    Notes:
358c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
359c8228a4eSBarry Smith 
360c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
36108405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
36208405cd6SLois Curfman McInnes .vb
36308405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
36408405cd6SLois Curfman McInnes       if (!(it % 2)) {
36508405cd6SLois Curfman McInnes         [compute Jacobian here]
36608405cd6SLois Curfman McInnes       }
36708405cd6SLois Curfman McInnes .ve
368c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
36908405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
370c8228a4eSBarry Smith 
37136851e7fSLois Curfman McInnes    Level: intermediate
37236851e7fSLois Curfman McInnes 
3739b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3749b94acceSBarry Smith @*/
37563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter)
3769b94acceSBarry Smith {
3773a40ed3dSBarry Smith   PetscFunctionBegin;
3784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3794482741eSBarry Smith   PetscValidIntPointer(iter,2);
3809b94acceSBarry Smith   *iter = snes->iter;
3813a40ed3dSBarry Smith   PetscFunctionReturn(0);
3829b94acceSBarry Smith }
38374679c65SBarry Smith 
3844a2ae208SSatish Balay #undef __FUNCT__
3854a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
3869b94acceSBarry Smith /*@
3879b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3889b94acceSBarry Smith    with SNESSSetFunction().
3899b94acceSBarry Smith 
390c7afd0dbSLois Curfman McInnes    Collective on SNES
391c7afd0dbSLois Curfman McInnes 
3929b94acceSBarry Smith    Input Parameter:
3939b94acceSBarry Smith .  snes - SNES context
3949b94acceSBarry Smith 
3959b94acceSBarry Smith    Output Parameter:
3969b94acceSBarry Smith .  fnorm - 2-norm of function
3979b94acceSBarry Smith 
39836851e7fSLois Curfman McInnes    Level: intermediate
39936851e7fSLois Curfman McInnes 
4009b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
401a86d99e1SLois Curfman McInnes 
40208405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
4039b94acceSBarry Smith @*/
40463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
4059b94acceSBarry Smith {
4063a40ed3dSBarry Smith   PetscFunctionBegin;
4074482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4084482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
4099b94acceSBarry Smith   *fnorm = snes->norm;
4103a40ed3dSBarry Smith   PetscFunctionReturn(0);
4119b94acceSBarry Smith }
41274679c65SBarry Smith 
4134a2ae208SSatish Balay #undef __FUNCT__
4144a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
4159b94acceSBarry Smith /*@
4169b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4179b94acceSBarry Smith    attempted by the nonlinear solver.
4189b94acceSBarry Smith 
419c7afd0dbSLois Curfman McInnes    Not Collective
420c7afd0dbSLois Curfman McInnes 
4219b94acceSBarry Smith    Input Parameter:
4229b94acceSBarry Smith .  snes - SNES context
4239b94acceSBarry Smith 
4249b94acceSBarry Smith    Output Parameter:
4259b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4269b94acceSBarry Smith 
427c96a6f78SLois Curfman McInnes    Notes:
428c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
429c96a6f78SLois Curfman McInnes 
43036851e7fSLois Curfman McInnes    Level: intermediate
43136851e7fSLois Curfman McInnes 
4329b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4339b94acceSBarry Smith @*/
43463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* nfails)
4359b94acceSBarry Smith {
4363a40ed3dSBarry Smith   PetscFunctionBegin;
4374482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4384482741eSBarry Smith   PetscValidIntPointer(nfails,2);
43950ffb88aSMatthew Knepley   *nfails = snes->numFailures;
44050ffb88aSMatthew Knepley   PetscFunctionReturn(0);
44150ffb88aSMatthew Knepley }
44250ffb88aSMatthew Knepley 
44350ffb88aSMatthew Knepley #undef __FUNCT__
44450ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
44550ffb88aSMatthew Knepley /*@
44650ffb88aSMatthew Knepley    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
44750ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
44850ffb88aSMatthew Knepley 
44950ffb88aSMatthew Knepley    Not Collective
45050ffb88aSMatthew Knepley 
45150ffb88aSMatthew Knepley    Input Parameters:
45250ffb88aSMatthew Knepley +  snes     - SNES context
45350ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
45450ffb88aSMatthew Knepley 
45550ffb88aSMatthew Knepley    Level: intermediate
45650ffb88aSMatthew Knepley 
45750ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
45850ffb88aSMatthew Knepley @*/
45963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt maxFails)
46050ffb88aSMatthew Knepley {
46150ffb88aSMatthew Knepley   PetscFunctionBegin;
4624482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
46350ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
46450ffb88aSMatthew Knepley   PetscFunctionReturn(0);
46550ffb88aSMatthew Knepley }
46650ffb88aSMatthew Knepley 
46750ffb88aSMatthew Knepley #undef __FUNCT__
46850ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
46950ffb88aSMatthew Knepley /*@
47050ffb88aSMatthew Knepley    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
47150ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
47250ffb88aSMatthew Knepley 
47350ffb88aSMatthew Knepley    Not Collective
47450ffb88aSMatthew Knepley 
47550ffb88aSMatthew Knepley    Input Parameter:
47650ffb88aSMatthew Knepley .  snes     - SNES context
47750ffb88aSMatthew Knepley 
47850ffb88aSMatthew Knepley    Output Parameter:
47950ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
48050ffb88aSMatthew Knepley 
48150ffb88aSMatthew Knepley    Level: intermediate
48250ffb88aSMatthew Knepley 
48350ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
48450ffb88aSMatthew Knepley @*/
48563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *maxFails)
48650ffb88aSMatthew Knepley {
48750ffb88aSMatthew Knepley   PetscFunctionBegin;
4884482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4894482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
49050ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
4913a40ed3dSBarry Smith   PetscFunctionReturn(0);
4929b94acceSBarry Smith }
493a847f771SSatish Balay 
4944a2ae208SSatish Balay #undef __FUNCT__
4954a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations"
496c96a6f78SLois Curfman McInnes /*@
497c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
498c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
499c96a6f78SLois Curfman McInnes 
500c7afd0dbSLois Curfman McInnes    Not Collective
501c7afd0dbSLois Curfman McInnes 
502c96a6f78SLois Curfman McInnes    Input Parameter:
503c96a6f78SLois Curfman McInnes .  snes - SNES context
504c96a6f78SLois Curfman McInnes 
505c96a6f78SLois Curfman McInnes    Output Parameter:
506c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
507c96a6f78SLois Curfman McInnes 
508c96a6f78SLois Curfman McInnes    Notes:
509c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
510c96a6f78SLois Curfman McInnes 
51136851e7fSLois Curfman McInnes    Level: intermediate
51236851e7fSLois Curfman McInnes 
513c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
514c96a6f78SLois Curfman McInnes @*/
51563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberLinearIterations(SNES snes,PetscInt* lits)
516c96a6f78SLois Curfman McInnes {
5173a40ed3dSBarry Smith   PetscFunctionBegin;
5184482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5194482741eSBarry Smith   PetscValidIntPointer(lits,2);
520c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5213a40ed3dSBarry Smith   PetscFunctionReturn(0);
522c96a6f78SLois Curfman McInnes }
523c96a6f78SLois Curfman McInnes 
5244a2ae208SSatish Balay #undef __FUNCT__
52594b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
52652baeb72SSatish Balay /*@
52794b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
5289b94acceSBarry Smith 
52994b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
530c7afd0dbSLois Curfman McInnes 
5319b94acceSBarry Smith    Input Parameter:
5329b94acceSBarry Smith .  snes - the SNES context
5339b94acceSBarry Smith 
5349b94acceSBarry Smith    Output Parameter:
53594b7f48cSBarry Smith .  ksp - the KSP context
5369b94acceSBarry Smith 
5379b94acceSBarry Smith    Notes:
53894b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
5399b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5409b94acceSBarry Smith    KSP and PC contexts as well.
5419b94acceSBarry Smith 
54236851e7fSLois Curfman McInnes    Level: beginner
54336851e7fSLois Curfman McInnes 
54494b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
5459b94acceSBarry Smith 
54694b7f48cSBarry Smith .seealso: KSPGetPC()
5479b94acceSBarry Smith @*/
54863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp)
5499b94acceSBarry Smith {
5503a40ed3dSBarry Smith   PetscFunctionBegin;
5514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5524482741eSBarry Smith   PetscValidPointer(ksp,2);
55394b7f48cSBarry Smith   *ksp = snes->ksp;
5543a40ed3dSBarry Smith   PetscFunctionReturn(0);
5559b94acceSBarry Smith }
55682bf6240SBarry Smith 
5574a2ae208SSatish Balay #undef __FUNCT__
5584a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
5596849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
560e24b481bSBarry Smith {
561e24b481bSBarry Smith   PetscFunctionBegin;
562e24b481bSBarry Smith   PetscFunctionReturn(0);
563e24b481bSBarry Smith }
564e24b481bSBarry Smith 
5659b94acceSBarry Smith /* -----------------------------------------------------------*/
5664a2ae208SSatish Balay #undef __FUNCT__
5674a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
56852baeb72SSatish Balay /*@
5699b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5709b94acceSBarry Smith 
571c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
572c7afd0dbSLois Curfman McInnes 
573c7afd0dbSLois Curfman McInnes    Input Parameters:
574c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
5759b94acceSBarry Smith 
5769b94acceSBarry Smith    Output Parameter:
5779b94acceSBarry Smith .  outsnes - the new SNES context
5789b94acceSBarry Smith 
579c7afd0dbSLois Curfman McInnes    Options Database Keys:
580c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
581c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
582c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
583c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
584c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
585c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
586c1f60f51SBarry Smith 
58736851e7fSLois Curfman McInnes    Level: beginner
58836851e7fSLois Curfman McInnes 
5899b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5909b94acceSBarry Smith 
5914b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES
5929b94acceSBarry Smith @*/
59363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes)
5949b94acceSBarry Smith {
595dfbe8321SBarry Smith   PetscErrorCode      ierr;
5969b94acceSBarry Smith   SNES                snes;
5979b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
59837fcc0dbSBarry Smith 
5993a40ed3dSBarry Smith   PetscFunctionBegin;
600ed1caa07SMatthew Knepley   PetscValidPointer(outsnes,2);
6018ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
6028ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
6038ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
6048ba1e511SMatthew Knepley #endif
6058ba1e511SMatthew Knepley 
60652e6d16bSBarry Smith   ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
607e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
6089b94acceSBarry Smith   snes->max_its           = 50;
6099750a799SBarry Smith   snes->max_funcs	  = 10000;
6109b94acceSBarry Smith   snes->norm		  = 0.0;
611b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
612b4874afaSBarry Smith   snes->ttol              = 0.0;
61370441072SBarry Smith   snes->abstol		  = 1.e-50;
6149b94acceSBarry Smith   snes->xtol		  = 1.e-8;
6154b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
6169b94acceSBarry Smith   snes->nfuncs            = 0;
61750ffb88aSMatthew Knepley   snes->numFailures       = 0;
61850ffb88aSMatthew Knepley   snes->maxFailures       = 1;
6197a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
620639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6219b94acceSBarry Smith   snes->data              = 0;
6229b94acceSBarry Smith   snes->view              = 0;
6234dc4c822SBarry Smith   snes->setupcalled       = PETSC_FALSE;
624186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
6256f24a144SLois Curfman McInnes   snes->vwork             = 0;
6266f24a144SLois Curfman McInnes   snes->nwork             = 0;
627758f92a0SBarry Smith   snes->conv_hist_len     = 0;
628758f92a0SBarry Smith   snes->conv_hist_max     = 0;
629758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
630758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
631758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
632184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
6339b94acceSBarry Smith 
6349b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
635b0a32e0cSBarry Smith   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
63652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));CHKERRQ(ierr);
6379b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6389b94acceSBarry Smith   kctx->version     = 2;
6399b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6409b94acceSBarry Smith                              this was too large for some test cases */
6419b94acceSBarry Smith   kctx->rtol_last   = 0;
6429b94acceSBarry Smith   kctx->rtol_max    = .9;
6439b94acceSBarry Smith   kctx->gamma       = 1.0;
6449b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6459b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6469b94acceSBarry Smith   kctx->threshold   = .1;
6479b94acceSBarry Smith   kctx->lresid_last = 0;
6489b94acceSBarry Smith   kctx->norm_last   = 0;
6499b94acceSBarry Smith 
65094b7f48cSBarry Smith   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
65152e6d16bSBarry Smith   ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
6525334005bSBarry Smith 
6539b94acceSBarry Smith   *outsnes = snes;
65400036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
6553a40ed3dSBarry Smith   PetscFunctionReturn(0);
6569b94acceSBarry Smith }
6579b94acceSBarry Smith 
6584a2ae208SSatish Balay #undef __FUNCT__
6594a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
6609b94acceSBarry Smith /*@C
6619b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6629b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6639b94acceSBarry Smith    equations.
6649b94acceSBarry Smith 
665fee21e36SBarry Smith    Collective on SNES
666fee21e36SBarry Smith 
667c7afd0dbSLois Curfman McInnes    Input Parameters:
668c7afd0dbSLois Curfman McInnes +  snes - the SNES context
669c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
670c7afd0dbSLois Curfman McInnes .  r - vector to store function value
671c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
672c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6739b94acceSBarry Smith 
674c7afd0dbSLois Curfman McInnes    Calling sequence of func:
6758d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
676c7afd0dbSLois Curfman McInnes 
677313e4042SLois Curfman McInnes .  f - function vector
678c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6799b94acceSBarry Smith 
6809b94acceSBarry Smith    Notes:
6819b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6829b94acceSBarry Smith $      f'(x) x = -f(x),
683c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6849b94acceSBarry Smith 
68536851e7fSLois Curfman McInnes    Level: beginner
68636851e7fSLois Curfman McInnes 
6879b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6889b94acceSBarry Smith 
689a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6909b94acceSBarry Smith @*/
69163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
6929b94acceSBarry Smith {
6933a40ed3dSBarry Smith   PetscFunctionBegin;
6944482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6954482741eSBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
696c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,r,2);
697184914b5SBarry Smith 
6989b94acceSBarry Smith   snes->computefunction     = func;
6999b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7009b94acceSBarry Smith   snes->funP                = ctx;
7013a40ed3dSBarry Smith   PetscFunctionReturn(0);
7029b94acceSBarry Smith }
7039b94acceSBarry Smith 
7043ab0aad5SBarry Smith /* --------------------------------------------------------------- */
7053ab0aad5SBarry Smith #undef __FUNCT__
7063ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs"
7073ab0aad5SBarry Smith /*@C
7083ab0aad5SBarry Smith    SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set
7093ab0aad5SBarry Smith    it assumes a zero right hand side.
7103ab0aad5SBarry Smith 
7113ab0aad5SBarry Smith    Collective on SNES
7123ab0aad5SBarry Smith 
7133ab0aad5SBarry Smith    Input Parameters:
7143ab0aad5SBarry Smith +  snes - the SNES context
7153ab0aad5SBarry Smith -  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7163ab0aad5SBarry Smith 
7173ab0aad5SBarry Smith    Level: intermediate
7183ab0aad5SBarry Smith 
7193ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side
7203ab0aad5SBarry Smith 
7211096aae1SMatthew Knepley .seealso: SNESGetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7223ab0aad5SBarry Smith @*/
72363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetRhs(SNES snes,Vec rhs)
7243ab0aad5SBarry Smith {
725dfbe8321SBarry Smith   PetscErrorCode ierr;
7263ab0aad5SBarry Smith 
7273ab0aad5SBarry Smith   PetscFunctionBegin;
7283ab0aad5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7293ab0aad5SBarry Smith   if (rhs) {
7303ab0aad5SBarry Smith     PetscValidHeaderSpecific(rhs,VEC_COOKIE,2);
7313ab0aad5SBarry Smith     PetscCheckSameComm(snes,1,rhs,2);
7323ab0aad5SBarry Smith     ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr);
7333ab0aad5SBarry Smith   }
7343ab0aad5SBarry Smith   if (snes->afine) {
7353ab0aad5SBarry Smith     ierr = VecDestroy(snes->afine);CHKERRQ(ierr);
7363ab0aad5SBarry Smith   }
7373ab0aad5SBarry Smith   snes->afine = rhs;
7383ab0aad5SBarry Smith   PetscFunctionReturn(0);
7393ab0aad5SBarry Smith }
7403ab0aad5SBarry Smith 
7414a2ae208SSatish Balay #undef __FUNCT__
7421096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs"
7431096aae1SMatthew Knepley /*@C
7441096aae1SMatthew Knepley    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
7451096aae1SMatthew Knepley    it assumes a zero right hand side.
7461096aae1SMatthew Knepley 
7471096aae1SMatthew Knepley    Collective on SNES
7481096aae1SMatthew Knepley 
7491096aae1SMatthew Knepley    Input Parameter:
7501096aae1SMatthew Knepley .  snes - the SNES context
7511096aae1SMatthew Knepley 
7521096aae1SMatthew Knepley    Output Parameter:
7531096aae1SMatthew Knepley .  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7541096aae1SMatthew Knepley 
7551096aae1SMatthew Knepley    Level: intermediate
7561096aae1SMatthew Knepley 
7571096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side
7581096aae1SMatthew Knepley 
7591096aae1SMatthew Knepley .seealso: SNESSetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7601096aae1SMatthew Knepley @*/
7611096aae1SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs)
7621096aae1SMatthew Knepley {
7631096aae1SMatthew Knepley   PetscFunctionBegin;
7641096aae1SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7651096aae1SMatthew Knepley   PetscValidPointer(rhs,2);
7661096aae1SMatthew Knepley   *rhs = snes->afine;
7671096aae1SMatthew Knepley   PetscFunctionReturn(0);
7681096aae1SMatthew Knepley }
7691096aae1SMatthew Knepley 
7701096aae1SMatthew Knepley #undef __FUNCT__
7714a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
7729b94acceSBarry Smith /*@
77336851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
7749b94acceSBarry Smith                          SNESSetFunction().
7759b94acceSBarry Smith 
776c7afd0dbSLois Curfman McInnes    Collective on SNES
777c7afd0dbSLois Curfman McInnes 
7789b94acceSBarry Smith    Input Parameters:
779c7afd0dbSLois Curfman McInnes +  snes - the SNES context
780c7afd0dbSLois Curfman McInnes -  x - input vector
7819b94acceSBarry Smith 
7829b94acceSBarry Smith    Output Parameter:
7833638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7849b94acceSBarry Smith 
7851bffabb2SLois Curfman McInnes    Notes:
78636851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
78736851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
78836851e7fSLois Curfman McInnes    themselves.
78936851e7fSLois Curfman McInnes 
79036851e7fSLois Curfman McInnes    Level: developer
79136851e7fSLois Curfman McInnes 
7929b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7939b94acceSBarry Smith 
794a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7959b94acceSBarry Smith @*/
79663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y)
7979b94acceSBarry Smith {
798dfbe8321SBarry Smith   PetscErrorCode ierr;
7999b94acceSBarry Smith 
8003a40ed3dSBarry Smith   PetscFunctionBegin;
8014482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8024482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
8034482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
804c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
805c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
806184914b5SBarry Smith 
807d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8081096aae1SMatthew Knepley   if (snes->computefunction) {
809d64ed03dSBarry Smith     PetscStackPush("SNES user function");
810e9a2bbcdSBarry Smith     CHKMEMQ;
81119717074SBarry Smith     ierr = (*snes->computefunction)(snes,x,y,snes->funP);
812e9a2bbcdSBarry Smith     CHKMEMQ;
813d64ed03dSBarry Smith     PetscStackPop;
814d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
81519717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
81619717074SBarry Smith     }
817d5e45103SBarry Smith     CHKERRQ(ierr);
8181096aae1SMatthew Knepley   } else if (snes->afine) {
8191096aae1SMatthew Knepley     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
8201096aae1SMatthew Knepley   } else {
8211096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
8221096aae1SMatthew Knepley   }
8233ab0aad5SBarry Smith   if (snes->afine) {
824016dedfbSBarry Smith     ierr = VecAXPY(y,-1.0,snes->afine);CHKERRQ(ierr);
8253ab0aad5SBarry Smith   }
826ae3c334cSLois Curfman McInnes   snes->nfuncs++;
827d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8283a40ed3dSBarry Smith   PetscFunctionReturn(0);
8299b94acceSBarry Smith }
8309b94acceSBarry Smith 
8314a2ae208SSatish Balay #undef __FUNCT__
8324a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
83362fef451SLois Curfman McInnes /*@
83462fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
83562fef451SLois Curfman McInnes    set with SNESSetJacobian().
83662fef451SLois Curfman McInnes 
837c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
838c7afd0dbSLois Curfman McInnes 
83962fef451SLois Curfman McInnes    Input Parameters:
840c7afd0dbSLois Curfman McInnes +  snes - the SNES context
841c7afd0dbSLois Curfman McInnes -  x - input vector
84262fef451SLois Curfman McInnes 
84362fef451SLois Curfman McInnes    Output Parameters:
844c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
84562fef451SLois Curfman McInnes .  B - optional preconditioning matrix
846c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
847fee21e36SBarry Smith 
84862fef451SLois Curfman McInnes    Notes:
84962fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
85062fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
85162fef451SLois Curfman McInnes 
85294b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
853dc5a77f8SLois Curfman McInnes    flag parameter.
85462fef451SLois Curfman McInnes 
85536851e7fSLois Curfman McInnes    Level: developer
85636851e7fSLois Curfman McInnes 
85762fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
85862fef451SLois Curfman McInnes 
85994b7f48cSBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators()
86062fef451SLois Curfman McInnes @*/
86163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8629b94acceSBarry Smith {
863dfbe8321SBarry Smith   PetscErrorCode ierr;
8643a40ed3dSBarry Smith 
8653a40ed3dSBarry Smith   PetscFunctionBegin;
8664482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8674482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
8684482741eSBarry Smith   PetscValidPointer(flg,5);
869c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
8703a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
871d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
872c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
873d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
874dc67937bSBarry Smith   CHKMEMQ;
8759b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
876dc67937bSBarry Smith   CHKMEMQ;
877d64ed03dSBarry Smith   PetscStackPop;
878d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
8796d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
8804482741eSBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
8814482741eSBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE,4);
8823a40ed3dSBarry Smith   PetscFunctionReturn(0);
8839b94acceSBarry Smith }
8849b94acceSBarry Smith 
8854a2ae208SSatish Balay #undef __FUNCT__
8864a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
8879b94acceSBarry Smith /*@C
8889b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
889044dda88SLois Curfman McInnes    location to store the matrix.
8909b94acceSBarry Smith 
891c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
892c7afd0dbSLois Curfman McInnes 
8939b94acceSBarry Smith    Input Parameters:
894c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8959b94acceSBarry Smith .  A - Jacobian matrix
8969b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8979b94acceSBarry Smith .  func - Jacobian evaluation routine
898c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
8992cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
9009b94acceSBarry Smith 
9019b94acceSBarry Smith    Calling sequence of func:
9028d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
9039b94acceSBarry Smith 
904c7afd0dbSLois Curfman McInnes +  x - input vector
9059b94acceSBarry Smith .  A - Jacobian matrix
9069b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
907ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
90894b7f48cSBarry Smith    structure (same as flag in KSPSetOperators())
909c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
9109b94acceSBarry Smith 
9119b94acceSBarry Smith    Notes:
91294b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
9132cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
914ac21db08SLois Curfman McInnes 
915ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9169b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9179b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9189b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9199b94acceSBarry Smith    throughout the global iterations.
9209b94acceSBarry Smith 
92136851e7fSLois Curfman McInnes    Level: beginner
92236851e7fSLois Curfman McInnes 
9239b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9249b94acceSBarry Smith 
9253960baedSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor()
9269b94acceSBarry Smith @*/
92763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
9289b94acceSBarry Smith {
929dfbe8321SBarry Smith   PetscErrorCode ierr;
9303a7fca6bSBarry Smith 
9313a40ed3dSBarry Smith   PetscFunctionBegin;
9324482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9334482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
9344482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
935c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
936c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
9373a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
9383a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
9393a7fca6bSBarry Smith   if (A) {
9403a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
9419b94acceSBarry Smith     snes->jacobian = A;
9423a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
9433a7fca6bSBarry Smith   }
9443a7fca6bSBarry Smith   if (B) {
9453a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
9469b94acceSBarry Smith     snes->jacobian_pre = B;
9473a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
9483a7fca6bSBarry Smith   }
94969a612a9SBarry Smith   ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9503a40ed3dSBarry Smith   PetscFunctionReturn(0);
9519b94acceSBarry Smith }
95262fef451SLois Curfman McInnes 
9534a2ae208SSatish Balay #undef __FUNCT__
9544a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
955c2aafc4cSSatish Balay /*@C
956b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
957b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
958b4fd4287SBarry Smith 
959c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
960c7afd0dbSLois Curfman McInnes 
961b4fd4287SBarry Smith    Input Parameter:
962b4fd4287SBarry Smith .  snes - the nonlinear solver context
963b4fd4287SBarry Smith 
964b4fd4287SBarry Smith    Output Parameters:
965c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
966b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
96770e92668SMatthew Knepley .  func - location to put Jacobian function (or PETSC_NULL)
96870e92668SMatthew Knepley -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
969fee21e36SBarry Smith 
97036851e7fSLois Curfman McInnes    Level: advanced
97136851e7fSLois Curfman McInnes 
972b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
973b4fd4287SBarry Smith @*/
97470e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
975b4fd4287SBarry Smith {
9763a40ed3dSBarry Smith   PetscFunctionBegin;
9774482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
978b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
979b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
98000036973SBarry Smith   if (func) *func = snes->computejacobian;
98170e92668SMatthew Knepley   if (ctx)  *ctx  = snes->jacP;
9823a40ed3dSBarry Smith   PetscFunctionReturn(0);
983b4fd4287SBarry Smith }
984b4fd4287SBarry Smith 
9859b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
98663dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
9879b94acceSBarry Smith 
9884a2ae208SSatish Balay #undef __FUNCT__
9894a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
9909b94acceSBarry Smith /*@
9919b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
992272ac6f2SLois Curfman McInnes    of a nonlinear solver.
9939b94acceSBarry Smith 
994fee21e36SBarry Smith    Collective on SNES
995fee21e36SBarry Smith 
996c7afd0dbSLois Curfman McInnes    Input Parameters:
99770e92668SMatthew Knepley .  snes - the SNES context
998c7afd0dbSLois Curfman McInnes 
999272ac6f2SLois Curfman McInnes    Notes:
1000272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1001272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1002272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1003272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1004272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1005272ac6f2SLois Curfman McInnes 
100636851e7fSLois Curfman McInnes    Level: advanced
100736851e7fSLois Curfman McInnes 
10089b94acceSBarry Smith .keywords: SNES, nonlinear, setup
10099b94acceSBarry Smith 
10109b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
10119b94acceSBarry Smith @*/
101270e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
10139b94acceSBarry Smith {
1014dfbe8321SBarry Smith   PetscErrorCode ierr;
10154b27c08aSLois Curfman McInnes   PetscTruth     flg, iseqtr;
10163a40ed3dSBarry Smith 
10173a40ed3dSBarry Smith   PetscFunctionBegin;
10184482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10194dc4c822SBarry Smith   if (snes->setupcalled) PetscFunctionReturn(0);
10209b94acceSBarry Smith 
1021b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1022c1f60f51SBarry Smith   /*
1023c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1024dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1025c1f60f51SBarry Smith   */
1026c1f60f51SBarry Smith   if (flg) {
1027c1f60f51SBarry Smith     Mat J;
10285a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10295a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
103063ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr);
10313a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
10323a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1033c1f60f51SBarry Smith   }
103445fc7adcSBarry Smith 
1035c8d321e0SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE)
103645fc7adcSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
103745fc7adcSBarry Smith   if (flg) {
103845fc7adcSBarry Smith     Mat J;
103945fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
104045fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
104145fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
104245fc7adcSBarry Smith   }
104332a4b47aSMatthew Knepley #endif
104445fc7adcSBarry Smith 
1045b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1046c1f60f51SBarry Smith   /*
1047dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1048c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1049c1f60f51SBarry Smith    */
105031615d04SLois Curfman McInnes   if (flg) {
1051272ac6f2SLois Curfman McInnes     Mat  J;
1052b5d62d44SBarry Smith     KSP ksp;
105394b7f48cSBarry Smith     PC   pc;
1054f3ef73ceSBarry Smith 
10555a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10563a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
105763ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr);
10583a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
10593a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
10603a7fca6bSBarry Smith 
1061f3ef73ceSBarry Smith     /* force no preconditioner */
106294b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1063b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1064a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1065a9815358SBarry Smith     if (!flg) {
1066f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1067272ac6f2SLois Curfman McInnes     }
1068a9815358SBarry Smith   }
1069f3ef73ceSBarry Smith 
10701096aae1SMatthew Knepley   if (!snes->vec_func && !snes->afine) {
10711096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10721096aae1SMatthew Knepley   }
10731096aae1SMatthew Knepley   if (!snes->computefunction && !snes->afine) {
10741096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10751096aae1SMatthew Knepley   }
107629bbc08cSBarry Smith   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1077a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
107829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1079a8c6a408SBarry Smith   }
1080a703fe33SLois Curfman McInnes 
1081a703fe33SLois Curfman McInnes   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
10824b27c08aSLois Curfman McInnes   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
10836831982aSBarry Smith   if (snes->ksp_ewconv && !iseqtr) {
108494b7f48cSBarry Smith     KSP ksp;
108594b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1086186905e3SBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1087a703fe33SLois Curfman McInnes   }
10884b27c08aSLois Curfman McInnes 
1089a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
10907aaed0d8SBarry Smith   snes->setupcalled = PETSC_TRUE;
10913a40ed3dSBarry Smith   PetscFunctionReturn(0);
10929b94acceSBarry Smith }
10939b94acceSBarry Smith 
10944a2ae208SSatish Balay #undef __FUNCT__
10954a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
109652baeb72SSatish Balay /*@
10979b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
10989b94acceSBarry Smith    with SNESCreate().
10999b94acceSBarry Smith 
1100c7afd0dbSLois Curfman McInnes    Collective on SNES
1101c7afd0dbSLois Curfman McInnes 
11029b94acceSBarry Smith    Input Parameter:
11039b94acceSBarry Smith .  snes - the SNES context
11049b94acceSBarry Smith 
110536851e7fSLois Curfman McInnes    Level: beginner
110636851e7fSLois Curfman McInnes 
11079b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
11089b94acceSBarry Smith 
110963a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
11109b94acceSBarry Smith @*/
111163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
11129b94acceSBarry Smith {
11136849ba73SBarry Smith   PetscErrorCode ierr;
11143a40ed3dSBarry Smith 
11153a40ed3dSBarry Smith   PetscFunctionBegin;
11164482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
11173a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1118d4bb536fSBarry Smith 
1119be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
11200f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1121be0abb6dSBarry Smith 
1122e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1123606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
11243a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
11253a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
11263ab0aad5SBarry Smith   if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);}
112794b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1128522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1129d952e501SBarry Smith   ierr = SNESClearMonitor(snes);CHKERRQ(ierr);
1130a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
11313a40ed3dSBarry Smith   PetscFunctionReturn(0);
11329b94acceSBarry Smith }
11339b94acceSBarry Smith 
11349b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11359b94acceSBarry Smith 
11364a2ae208SSatish Balay #undef __FUNCT__
11374a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
11389b94acceSBarry Smith /*@
1139d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11409b94acceSBarry Smith 
1141c7afd0dbSLois Curfman McInnes    Collective on SNES
1142c7afd0dbSLois Curfman McInnes 
11439b94acceSBarry Smith    Input Parameters:
1144c7afd0dbSLois Curfman McInnes +  snes - the SNES context
114570441072SBarry Smith .  abstol - absolute convergence tolerance
114633174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
114733174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
114833174efeSLois Curfman McInnes            of the change in the solution between steps
114933174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1150c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1151fee21e36SBarry Smith 
115233174efeSLois Curfman McInnes    Options Database Keys:
115370441072SBarry Smith +    -snes_atol <abstol> - Sets abstol
1154c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1155c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1156c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1157c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
11589b94acceSBarry Smith 
1159d7a720efSLois Curfman McInnes    Notes:
11609b94acceSBarry Smith    The default maximum number of iterations is 50.
11619b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11629b94acceSBarry Smith 
116336851e7fSLois Curfman McInnes    Level: intermediate
116436851e7fSLois Curfman McInnes 
116533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11669b94acceSBarry Smith 
11672492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance()
11689b94acceSBarry Smith @*/
116963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
11709b94acceSBarry Smith {
11713a40ed3dSBarry Smith   PetscFunctionBegin;
11724482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
117370441072SBarry Smith   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1174d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1175d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1176d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1177d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11783a40ed3dSBarry Smith   PetscFunctionReturn(0);
11799b94acceSBarry Smith }
11809b94acceSBarry Smith 
11814a2ae208SSatish Balay #undef __FUNCT__
11824a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
11839b94acceSBarry Smith /*@
118433174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
118533174efeSLois Curfman McInnes 
1186c7afd0dbSLois Curfman McInnes    Not Collective
1187c7afd0dbSLois Curfman McInnes 
118833174efeSLois Curfman McInnes    Input Parameters:
1189c7afd0dbSLois Curfman McInnes +  snes - the SNES context
119070441072SBarry Smith .  abstol - absolute convergence tolerance
119133174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
119233174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
119333174efeSLois Curfman McInnes            of the change in the solution between steps
119433174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1195c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1196fee21e36SBarry Smith 
119733174efeSLois Curfman McInnes    Notes:
119833174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
119933174efeSLois Curfman McInnes 
120036851e7fSLois Curfman McInnes    Level: intermediate
120136851e7fSLois Curfman McInnes 
120233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
120333174efeSLois Curfman McInnes 
120433174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
120533174efeSLois Curfman McInnes @*/
120663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
120733174efeSLois Curfman McInnes {
12083a40ed3dSBarry Smith   PetscFunctionBegin;
12094482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
121070441072SBarry Smith   if (abstol)  *abstol  = snes->abstol;
121133174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
121233174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
121333174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
121433174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
12153a40ed3dSBarry Smith   PetscFunctionReturn(0);
121633174efeSLois Curfman McInnes }
121733174efeSLois Curfman McInnes 
12184a2ae208SSatish Balay #undef __FUNCT__
12194a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
122033174efeSLois Curfman McInnes /*@
12219b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12229b94acceSBarry Smith 
1223fee21e36SBarry Smith    Collective on SNES
1224fee21e36SBarry Smith 
1225c7afd0dbSLois Curfman McInnes    Input Parameters:
1226c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1227c7afd0dbSLois Curfman McInnes -  tol - tolerance
1228c7afd0dbSLois Curfman McInnes 
12299b94acceSBarry Smith    Options Database Key:
1230c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
12319b94acceSBarry Smith 
123236851e7fSLois Curfman McInnes    Level: intermediate
123336851e7fSLois Curfman McInnes 
12349b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12359b94acceSBarry Smith 
12362492ecdbSBarry Smith .seealso: SNESSetTolerances()
12379b94acceSBarry Smith @*/
123863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
12399b94acceSBarry Smith {
12403a40ed3dSBarry Smith   PetscFunctionBegin;
12414482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12429b94acceSBarry Smith   snes->deltatol = tol;
12433a40ed3dSBarry Smith   PetscFunctionReturn(0);
12449b94acceSBarry Smith }
12459b94acceSBarry Smith 
1246df9fa365SBarry Smith /*
1247df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1248df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1249df9fa365SBarry Smith    macros instead of functions
1250df9fa365SBarry Smith */
12514a2ae208SSatish Balay #undef __FUNCT__
12524a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
125363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1254ce1608b8SBarry Smith {
1255dfbe8321SBarry Smith   PetscErrorCode ierr;
1256ce1608b8SBarry Smith 
1257ce1608b8SBarry Smith   PetscFunctionBegin;
12584482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1259ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1260ce1608b8SBarry Smith   PetscFunctionReturn(0);
1261ce1608b8SBarry Smith }
1262ce1608b8SBarry Smith 
12634a2ae208SSatish Balay #undef __FUNCT__
12644a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
126563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1266df9fa365SBarry Smith {
1267dfbe8321SBarry Smith   PetscErrorCode ierr;
1268df9fa365SBarry Smith 
1269df9fa365SBarry Smith   PetscFunctionBegin;
1270df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1271df9fa365SBarry Smith   PetscFunctionReturn(0);
1272df9fa365SBarry Smith }
1273df9fa365SBarry Smith 
12744a2ae208SSatish Balay #undef __FUNCT__
12754a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
127663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw)
1277df9fa365SBarry Smith {
1278dfbe8321SBarry Smith   PetscErrorCode ierr;
1279df9fa365SBarry Smith 
1280df9fa365SBarry Smith   PetscFunctionBegin;
1281df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1282df9fa365SBarry Smith   PetscFunctionReturn(0);
1283df9fa365SBarry Smith }
1284df9fa365SBarry Smith 
12859b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12869b94acceSBarry Smith 
12874a2ae208SSatish Balay #undef __FUNCT__
12884a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
12899b94acceSBarry Smith /*@C
12905cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
12919b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12929b94acceSBarry Smith    progress.
12939b94acceSBarry Smith 
1294fee21e36SBarry Smith    Collective on SNES
1295fee21e36SBarry Smith 
1296c7afd0dbSLois Curfman McInnes    Input Parameters:
1297c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1298c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1299b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1300*e8105e01SRichard Katz           monitor routine (use PETSC_NULL if no context is desired)
1301b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1302b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
13039b94acceSBarry Smith 
1304c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1305a7cc72afSBarry Smith $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1306c7afd0dbSLois Curfman McInnes 
1307c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1308c7afd0dbSLois Curfman McInnes .    its - iteration number
1309c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
131040a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
13119b94acceSBarry Smith 
13129665c990SLois Curfman McInnes    Options Database Keys:
1313c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1314c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1315c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1316c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1317c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1318c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1319c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1320c7afd0dbSLois Curfman McInnes                             the options database.
13219665c990SLois Curfman McInnes 
1322639f9d9dSBarry Smith    Notes:
13236bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13246bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13256bc08f3fSLois Curfman McInnes    order in which they were set.
1326639f9d9dSBarry Smith 
132736851e7fSLois Curfman McInnes    Level: intermediate
132836851e7fSLois Curfman McInnes 
13299b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13309b94acceSBarry Smith 
13315cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
13329b94acceSBarry Smith @*/
133363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
13349b94acceSBarry Smith {
13353a40ed3dSBarry Smith   PetscFunctionBegin;
13364482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1337639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
133829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1339639f9d9dSBarry Smith   }
1340639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1341b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1342639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13433a40ed3dSBarry Smith   PetscFunctionReturn(0);
13449b94acceSBarry Smith }
13459b94acceSBarry Smith 
13464a2ae208SSatish Balay #undef __FUNCT__
13474a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
13485cd90555SBarry Smith /*@C
13495cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
13505cd90555SBarry Smith 
1351c7afd0dbSLois Curfman McInnes    Collective on SNES
1352c7afd0dbSLois Curfman McInnes 
13535cd90555SBarry Smith    Input Parameters:
13545cd90555SBarry Smith .  snes - the SNES context
13555cd90555SBarry Smith 
13561a480d89SAdministrator    Options Database Key:
1357c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1358c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1359c7afd0dbSLois Curfman McInnes     set via the options database
13605cd90555SBarry Smith 
13615cd90555SBarry Smith    Notes:
13625cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
13635cd90555SBarry Smith 
136436851e7fSLois Curfman McInnes    Level: intermediate
136536851e7fSLois Curfman McInnes 
13665cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
13675cd90555SBarry Smith 
13685cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
13695cd90555SBarry Smith @*/
137063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes)
13715cd90555SBarry Smith {
1372d952e501SBarry Smith   PetscErrorCode ierr;
1373d952e501SBarry Smith   PetscInt       i;
1374d952e501SBarry Smith 
13755cd90555SBarry Smith   PetscFunctionBegin;
13764482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1377d952e501SBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1378d952e501SBarry Smith     if (snes->monitordestroy[i]) {
1379d952e501SBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1380d952e501SBarry Smith     }
1381d952e501SBarry Smith   }
13825cd90555SBarry Smith   snes->numbermonitors = 0;
13835cd90555SBarry Smith   PetscFunctionReturn(0);
13845cd90555SBarry Smith }
13855cd90555SBarry Smith 
13864a2ae208SSatish Balay #undef __FUNCT__
13874a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
13889b94acceSBarry Smith /*@C
13899b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13909b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13919b94acceSBarry Smith 
1392fee21e36SBarry Smith    Collective on SNES
1393fee21e36SBarry Smith 
1394c7afd0dbSLois Curfman McInnes    Input Parameters:
1395c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1396c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1397c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1398c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
13999b94acceSBarry Smith 
1400c7afd0dbSLois Curfman McInnes    Calling sequence of func:
14016849ba73SBarry Smith $     PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1402c7afd0dbSLois Curfman McInnes 
1403c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1404c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1405184914b5SBarry Smith .    reason - reason for convergence/divergence
1406c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
14074b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
14084b27c08aSLois Curfman McInnes -    f - 2-norm of function
14099b94acceSBarry Smith 
141036851e7fSLois Curfman McInnes    Level: advanced
141136851e7fSLois Curfman McInnes 
14129b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
14139b94acceSBarry Smith 
14144b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR()
14159b94acceSBarry Smith @*/
141663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
14179b94acceSBarry Smith {
14183a40ed3dSBarry Smith   PetscFunctionBegin;
14194482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14209b94acceSBarry Smith   (snes)->converged = func;
14219b94acceSBarry Smith   (snes)->cnvP      = cctx;
14223a40ed3dSBarry Smith   PetscFunctionReturn(0);
14239b94acceSBarry Smith }
14249b94acceSBarry Smith 
14254a2ae208SSatish Balay #undef __FUNCT__
14264a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
142752baeb72SSatish Balay /*@
1428184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1429184914b5SBarry Smith 
1430184914b5SBarry Smith    Not Collective
1431184914b5SBarry Smith 
1432184914b5SBarry Smith    Input Parameter:
1433184914b5SBarry Smith .  snes - the SNES context
1434184914b5SBarry Smith 
1435184914b5SBarry Smith    Output Parameter:
1436e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1437184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1438184914b5SBarry Smith 
1439184914b5SBarry Smith    Level: intermediate
1440184914b5SBarry Smith 
1441184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1442184914b5SBarry Smith 
1443184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1444184914b5SBarry Smith 
14454b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1446184914b5SBarry Smith @*/
144763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1448184914b5SBarry Smith {
1449184914b5SBarry Smith   PetscFunctionBegin;
14504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14514482741eSBarry Smith   PetscValidPointer(reason,2);
1452184914b5SBarry Smith   *reason = snes->reason;
1453184914b5SBarry Smith   PetscFunctionReturn(0);
1454184914b5SBarry Smith }
1455184914b5SBarry Smith 
14564a2ae208SSatish Balay #undef __FUNCT__
14574a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1458c9005455SLois Curfman McInnes /*@
1459c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1460c9005455SLois Curfman McInnes 
1461fee21e36SBarry Smith    Collective on SNES
1462fee21e36SBarry Smith 
1463c7afd0dbSLois Curfman McInnes    Input Parameters:
1464c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1465c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1466cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1467758f92a0SBarry Smith .  na  - size of a and its
146864731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1469758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1470c7afd0dbSLois Curfman McInnes 
1471c9005455SLois Curfman McInnes    Notes:
14724b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1473c9005455SLois Curfman McInnes    at each step.
1474c9005455SLois Curfman McInnes 
1475c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1476c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1477c9005455SLois Curfman McInnes    during the section of code that is being timed.
1478c9005455SLois Curfman McInnes 
147936851e7fSLois Curfman McInnes    Level: intermediate
148036851e7fSLois Curfman McInnes 
1481c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1482758f92a0SBarry Smith 
148308405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1484758f92a0SBarry Smith 
1485c9005455SLois Curfman McInnes @*/
148663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset)
1487c9005455SLois Curfman McInnes {
14883a40ed3dSBarry Smith   PetscFunctionBegin;
14894482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14904482741eSBarry Smith   if (na) PetscValidScalarPointer(a,2);
1491c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1492758f92a0SBarry Smith   snes->conv_hist_its   = its;
1493758f92a0SBarry Smith   snes->conv_hist_max   = na;
1494758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1495758f92a0SBarry Smith   PetscFunctionReturn(0);
1496758f92a0SBarry Smith }
1497758f92a0SBarry Smith 
14984a2ae208SSatish Balay #undef __FUNCT__
14994a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
15000c4c9dddSBarry Smith /*@C
1501758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1502758f92a0SBarry Smith 
1503758f92a0SBarry Smith    Collective on SNES
1504758f92a0SBarry Smith 
1505758f92a0SBarry Smith    Input Parameter:
1506758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1507758f92a0SBarry Smith 
1508758f92a0SBarry Smith    Output Parameters:
1509758f92a0SBarry Smith .  a   - array to hold history
1510758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1511758f92a0SBarry Smith          negative if not converged) for each solve.
1512758f92a0SBarry Smith -  na  - size of a and its
1513758f92a0SBarry Smith 
1514758f92a0SBarry Smith    Notes:
1515758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1516758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1517758f92a0SBarry Smith 
1518758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1519758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1520758f92a0SBarry Smith    during the section of code that is being timed.
1521758f92a0SBarry Smith 
1522758f92a0SBarry Smith    Level: intermediate
1523758f92a0SBarry Smith 
1524758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1525758f92a0SBarry Smith 
1526758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1527758f92a0SBarry Smith 
1528758f92a0SBarry Smith @*/
152963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1530758f92a0SBarry Smith {
1531758f92a0SBarry Smith   PetscFunctionBegin;
15324482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1533758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1534758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1535758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
15363a40ed3dSBarry Smith   PetscFunctionReturn(0);
1537c9005455SLois Curfman McInnes }
1538c9005455SLois Curfman McInnes 
1539e74ef692SMatthew Knepley #undef __FUNCT__
1540e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1541ac226902SBarry Smith /*@C
154276b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
15437e4bb74cSBarry Smith   at the beginning o every iteration of the nonlinear solve. Specifically
15447e4bb74cSBarry Smith   it is called just before the Jacobian is "evaluated".
154576b2cf59SMatthew Knepley 
154676b2cf59SMatthew Knepley   Collective on SNES
154776b2cf59SMatthew Knepley 
154876b2cf59SMatthew Knepley   Input Parameters:
154976b2cf59SMatthew Knepley . snes - The nonlinear solver context
155076b2cf59SMatthew Knepley . func - The function
155176b2cf59SMatthew Knepley 
155276b2cf59SMatthew Knepley   Calling sequence of func:
1553b5d30489SBarry Smith . func (SNES snes, PetscInt step);
155476b2cf59SMatthew Knepley 
155576b2cf59SMatthew Knepley . step - The current step of the iteration
155676b2cf59SMatthew Knepley 
155776b2cf59SMatthew Knepley   Level: intermediate
155876b2cf59SMatthew Knepley 
155976b2cf59SMatthew Knepley .keywords: SNES, update
1560b5d30489SBarry Smith 
156176b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
156276b2cf59SMatthew Knepley @*/
156363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
156476b2cf59SMatthew Knepley {
156576b2cf59SMatthew Knepley   PetscFunctionBegin;
15664482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
156776b2cf59SMatthew Knepley   snes->update = func;
156876b2cf59SMatthew Knepley   PetscFunctionReturn(0);
156976b2cf59SMatthew Knepley }
157076b2cf59SMatthew Knepley 
1571e74ef692SMatthew Knepley #undef __FUNCT__
1572e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
157376b2cf59SMatthew Knepley /*@
157476b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
157576b2cf59SMatthew Knepley 
157676b2cf59SMatthew Knepley   Not collective
157776b2cf59SMatthew Knepley 
157876b2cf59SMatthew Knepley   Input Parameters:
157976b2cf59SMatthew Knepley . snes - The nonlinear solver context
158076b2cf59SMatthew Knepley . step - The current step of the iteration
158176b2cf59SMatthew Knepley 
1582205452f4SMatthew Knepley   Level: intermediate
1583205452f4SMatthew Knepley 
158476b2cf59SMatthew Knepley .keywords: SNES, update
158576b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
158676b2cf59SMatthew Knepley @*/
158763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
158876b2cf59SMatthew Knepley {
158976b2cf59SMatthew Knepley   PetscFunctionBegin;
159076b2cf59SMatthew Knepley   PetscFunctionReturn(0);
159176b2cf59SMatthew Knepley }
159276b2cf59SMatthew Knepley 
15934a2ae208SSatish Balay #undef __FUNCT__
15944a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
15959b94acceSBarry Smith /*
15969b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15979b94acceSBarry Smith    positive parameter delta.
15989b94acceSBarry Smith 
15999b94acceSBarry Smith     Input Parameters:
1600c7afd0dbSLois Curfman McInnes +   snes - the SNES context
16019b94acceSBarry Smith .   y - approximate solution of linear system
16029b94acceSBarry Smith .   fnorm - 2-norm of current function
1603c7afd0dbSLois Curfman McInnes -   delta - trust region size
16049b94acceSBarry Smith 
16059b94acceSBarry Smith     Output Parameters:
1606c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
16079b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
16089b94acceSBarry Smith     region, and exceeds zero otherwise.
1609c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
16109b94acceSBarry Smith 
16119b94acceSBarry Smith     Note:
16124b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
16139b94acceSBarry Smith     is set to be the maximum allowable step size.
16149b94acceSBarry Smith 
16159b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
16169b94acceSBarry Smith */
1617dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
16189b94acceSBarry Smith {
1619064f8208SBarry Smith   PetscReal      nrm;
1620ea709b57SSatish Balay   PetscScalar    cnorm;
1621dfbe8321SBarry Smith   PetscErrorCode ierr;
16223a40ed3dSBarry Smith 
16233a40ed3dSBarry Smith   PetscFunctionBegin;
16244482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16254482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1626c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1627184914b5SBarry Smith 
1628064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1629064f8208SBarry Smith   if (nrm > *delta) {
1630064f8208SBarry Smith      nrm = *delta/nrm;
1631064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1632064f8208SBarry Smith      cnorm = nrm;
16332dcb1b2aSMatthew Knepley      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
16349b94acceSBarry Smith      *ynorm = *delta;
16359b94acceSBarry Smith   } else {
16369b94acceSBarry Smith      *gpnorm = 0.0;
1637064f8208SBarry Smith      *ynorm = nrm;
16389b94acceSBarry Smith   }
16393a40ed3dSBarry Smith   PetscFunctionReturn(0);
16409b94acceSBarry Smith }
16419b94acceSBarry Smith 
16424a2ae208SSatish Balay #undef __FUNCT__
16434a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
16449b94acceSBarry Smith /*@
1645f69a0ea3SMatthew Knepley    SNESSolve - Solves a nonlinear system F(x) = b.
1646f69a0ea3SMatthew Knepley    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
16479b94acceSBarry Smith 
1648c7afd0dbSLois Curfman McInnes    Collective on SNES
1649c7afd0dbSLois Curfman McInnes 
1650b2002411SLois Curfman McInnes    Input Parameters:
1651c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1652f69a0ea3SMatthew Knepley .  b - the constant part of the equation, or PETSC_NULL to use zero.
165370e92668SMatthew Knepley -  x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution()
16549b94acceSBarry Smith 
1655b2002411SLois Curfman McInnes    Notes:
16568ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
16578ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16588ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16598ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16608ddd3da0SLois Curfman McInnes 
166136851e7fSLois Curfman McInnes    Level: beginner
166236851e7fSLois Curfman McInnes 
16639b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16649b94acceSBarry Smith 
166570e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution()
16669b94acceSBarry Smith @*/
1667f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
16689b94acceSBarry Smith {
1669dfbe8321SBarry Smith   PetscErrorCode ierr;
1670f1af5d2fSBarry Smith   PetscTruth     flg;
1671052efed2SBarry Smith 
16723a40ed3dSBarry Smith   PetscFunctionBegin;
16734482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16741302d50aSBarry Smith   if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
167574637425SBarry Smith 
1676f69a0ea3SMatthew Knepley   if (b) {
1677f69a0ea3SMatthew Knepley     ierr = SNESSetRhs(snes, b); CHKERRQ(ierr);
16781096aae1SMatthew Knepley     if (!snes->vec_func) {
16791096aae1SMatthew Knepley       Vec r;
16801096aae1SMatthew Knepley 
16811096aae1SMatthew Knepley       ierr = VecDuplicate(b, &r); CHKERRQ(ierr);
16821096aae1SMatthew Knepley       ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr);
16831096aae1SMatthew Knepley     }
1684f69a0ea3SMatthew Knepley   }
168570e92668SMatthew Knepley   if (x) {
1686f69a0ea3SMatthew Knepley     PetscValidHeaderSpecific(x,VEC_COOKIE,3);
1687f69a0ea3SMatthew Knepley     PetscCheckSameComm(snes,1,x,3);
168870e92668SMatthew Knepley   } else {
168970e92668SMatthew Knepley     ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr);
169070e92668SMatthew Knepley     if (!x) {
169170e92668SMatthew Knepley       ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr);
169270e92668SMatthew Knepley     }
169370e92668SMatthew Knepley   }
169470e92668SMatthew Knepley   snes->vec_sol = snes->vec_sol_always = x;
169570e92668SMatthew Knepley   if (!snes->setupcalled) {
169670e92668SMatthew Knepley     ierr = SNESSetUp(snes);CHKERRQ(ierr);
169770e92668SMatthew Knepley   }
1698abc0a331SBarry Smith   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
1699d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
170050ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1701d5e45103SBarry Smith 
1702d5e45103SBarry Smith   ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN);
1703d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
170438f152feSBarry Smith     /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */
170519717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr);
1706d5e45103SBarry Smith   } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) {
1707d5e45103SBarry Smith     /* translate exception into SNES not converged reason */
1708d5e45103SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
170938f152feSBarry Smith     ierr = 0;
171038f152feSBarry Smith   }
171138f152feSBarry Smith   CHKERRQ(ierr);
1712d5e45103SBarry Smith 
1713d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1714b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
17158b179ff8SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); }
1716da9b6338SBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1717da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
17185968eb51SBarry Smith   if (snes->printreason) {
17195968eb51SBarry Smith     if (snes->reason > 0) {
17209dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17215968eb51SBarry Smith     } else {
17229dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17235968eb51SBarry Smith     }
17245968eb51SBarry Smith   }
17255968eb51SBarry Smith 
17263a40ed3dSBarry Smith   PetscFunctionReturn(0);
17279b94acceSBarry Smith }
17289b94acceSBarry Smith 
17299b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
17309b94acceSBarry Smith 
17314a2ae208SSatish Balay #undef __FUNCT__
17324a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
173382bf6240SBarry Smith /*@C
17344b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
17359b94acceSBarry Smith 
1736fee21e36SBarry Smith    Collective on SNES
1737fee21e36SBarry Smith 
1738c7afd0dbSLois Curfman McInnes    Input Parameters:
1739c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1740454a90a3SBarry Smith -  type - a known method
1741c7afd0dbSLois Curfman McInnes 
1742c7afd0dbSLois Curfman McInnes    Options Database Key:
1743454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1744c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1745ae12b187SLois Curfman McInnes 
17469b94acceSBarry Smith    Notes:
1747e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
17484b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
1749c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17504b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
1751c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17529b94acceSBarry Smith 
1753ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1754ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1755ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1756ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1757ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1758ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1759ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1760ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1761ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1762b0a32e0cSBarry Smith   appropriate method.
176336851e7fSLois Curfman McInnes 
176436851e7fSLois Curfman McInnes   Level: intermediate
1765a703fe33SLois Curfman McInnes 
1766454a90a3SBarry Smith .keywords: SNES, set, type
1767435da068SBarry Smith 
1768435da068SBarry Smith .seealso: SNESType, SNESCreate()
1769435da068SBarry Smith 
17709b94acceSBarry Smith @*/
1771e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type)
17729b94acceSBarry Smith {
1773dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
17746831982aSBarry Smith   PetscTruth     match;
17753a40ed3dSBarry Smith 
17763a40ed3dSBarry Smith   PetscFunctionBegin;
17774482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17784482741eSBarry Smith   PetscValidCharPointer(type,2);
177982bf6240SBarry Smith 
17806831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
17810f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
178282bf6240SBarry Smith 
178382bf6240SBarry Smith   if (snes->setupcalled) {
17844dc4c822SBarry Smith     snes->setupcalled = PETSC_FALSE;
1785e1311b90SBarry Smith     ierr              = (*(snes)->destroy)(snes);CHKERRQ(ierr);
178682bf6240SBarry Smith     snes->data        = 0;
178782bf6240SBarry Smith   }
17887d1a2b2bSBarry Smith 
17899b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
179082bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1791b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1792958c9bccSBarry Smith   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
1793606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
179482bf6240SBarry Smith   snes->data = 0;
17953a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1796454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
17973a40ed3dSBarry Smith   PetscFunctionReturn(0);
17989b94acceSBarry Smith }
17999b94acceSBarry Smith 
1800a847f771SSatish Balay 
18019b94acceSBarry Smith /* --------------------------------------------------------------------- */
18024a2ae208SSatish Balay #undef __FUNCT__
18034a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
180452baeb72SSatish Balay /*@
18059b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1806f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
18079b94acceSBarry Smith 
1808fee21e36SBarry Smith    Not Collective
1809fee21e36SBarry Smith 
181036851e7fSLois Curfman McInnes    Level: advanced
181136851e7fSLois Curfman McInnes 
18129b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
18139b94acceSBarry Smith 
18149b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
18159b94acceSBarry Smith @*/
181663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
18179b94acceSBarry Smith {
1818dfbe8321SBarry Smith   PetscErrorCode ierr;
181982bf6240SBarry Smith 
18203a40ed3dSBarry Smith   PetscFunctionBegin;
182182bf6240SBarry Smith   if (SNESList) {
1822b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
182382bf6240SBarry Smith     SNESList = 0;
18249b94acceSBarry Smith   }
18254c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
18263a40ed3dSBarry Smith   PetscFunctionReturn(0);
18279b94acceSBarry Smith }
18289b94acceSBarry Smith 
18294a2ae208SSatish Balay #undef __FUNCT__
18304a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
18319b94acceSBarry Smith /*@C
18329a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
18339b94acceSBarry Smith 
1834c7afd0dbSLois Curfman McInnes    Not Collective
1835c7afd0dbSLois Curfman McInnes 
18369b94acceSBarry Smith    Input Parameter:
18374b0e389bSBarry Smith .  snes - nonlinear solver context
18389b94acceSBarry Smith 
18399b94acceSBarry Smith    Output Parameter:
18403a7fca6bSBarry Smith .  type - SNES method (a character string)
18419b94acceSBarry Smith 
184236851e7fSLois Curfman McInnes    Level: intermediate
184336851e7fSLois Curfman McInnes 
1844454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
18459b94acceSBarry Smith @*/
184663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type)
18479b94acceSBarry Smith {
18483a40ed3dSBarry Smith   PetscFunctionBegin;
18494482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18504482741eSBarry Smith   PetscValidPointer(type,2);
1851454a90a3SBarry Smith   *type = snes->type_name;
18523a40ed3dSBarry Smith   PetscFunctionReturn(0);
18539b94acceSBarry Smith }
18549b94acceSBarry Smith 
18554a2ae208SSatish Balay #undef __FUNCT__
18564a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
185752baeb72SSatish Balay /*@
18589b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
18599b94acceSBarry Smith    stored.
18609b94acceSBarry Smith 
1861c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1862c7afd0dbSLois Curfman McInnes 
18639b94acceSBarry Smith    Input Parameter:
18649b94acceSBarry Smith .  snes - the SNES context
18659b94acceSBarry Smith 
18669b94acceSBarry Smith    Output Parameter:
18679b94acceSBarry Smith .  x - the solution
18689b94acceSBarry Smith 
186970e92668SMatthew Knepley    Level: intermediate
187036851e7fSLois Curfman McInnes 
18719b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18729b94acceSBarry Smith 
187370e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
18749b94acceSBarry Smith @*/
187563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
18769b94acceSBarry Smith {
18773a40ed3dSBarry Smith   PetscFunctionBegin;
18784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18794482741eSBarry Smith   PetscValidPointer(x,2);
18809b94acceSBarry Smith   *x = snes->vec_sol_always;
18813a40ed3dSBarry Smith   PetscFunctionReturn(0);
18829b94acceSBarry Smith }
18839b94acceSBarry Smith 
18844a2ae208SSatish Balay #undef __FUNCT__
188570e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution"
18867e4bb74cSBarry Smith /*@
188770e92668SMatthew Knepley    SNESSetSolution - Sets the vector where the approximate solution is stored.
188870e92668SMatthew Knepley 
188970e92668SMatthew Knepley    Not Collective, but Vec is parallel if SNES is parallel
189070e92668SMatthew Knepley 
189170e92668SMatthew Knepley    Input Parameters:
189270e92668SMatthew Knepley +  snes - the SNES context
189370e92668SMatthew Knepley -  x - the solution
189470e92668SMatthew Knepley 
189570e92668SMatthew Knepley    Output Parameter:
189670e92668SMatthew Knepley 
189770e92668SMatthew Knepley    Level: intermediate
189870e92668SMatthew Knepley 
189942219521SBarry Smith    Notes: this is not normally used, rather one simply calls SNESSolve() with
190042219521SBarry Smith           the appropriate solution vector.
190142219521SBarry Smith 
190270e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution
190370e92668SMatthew Knepley 
190470e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
190570e92668SMatthew Knepley @*/
190670e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x)
190770e92668SMatthew Knepley {
190870e92668SMatthew Knepley   PetscFunctionBegin;
190970e92668SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
191070e92668SMatthew Knepley   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
191170e92668SMatthew Knepley   PetscCheckSameComm(snes,1,x,2);
191270e92668SMatthew Knepley   snes->vec_sol_always = x;
191370e92668SMatthew Knepley   PetscFunctionReturn(0);
191470e92668SMatthew Knepley }
191570e92668SMatthew Knepley 
191670e92668SMatthew Knepley #undef __FUNCT__
19174a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
191852baeb72SSatish Balay /*@
19199b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
19209b94acceSBarry Smith    stored.
19219b94acceSBarry Smith 
1922c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1923c7afd0dbSLois Curfman McInnes 
19249b94acceSBarry Smith    Input Parameter:
19259b94acceSBarry Smith .  snes - the SNES context
19269b94acceSBarry Smith 
19279b94acceSBarry Smith    Output Parameter:
19289b94acceSBarry Smith .  x - the solution update
19299b94acceSBarry Smith 
193036851e7fSLois Curfman McInnes    Level: advanced
193136851e7fSLois Curfman McInnes 
19329b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
19339b94acceSBarry Smith 
19349b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
19359b94acceSBarry Smith @*/
193663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
19379b94acceSBarry Smith {
19383a40ed3dSBarry Smith   PetscFunctionBegin;
19394482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19404482741eSBarry Smith   PetscValidPointer(x,2);
19419b94acceSBarry Smith   *x = snes->vec_sol_update_always;
19423a40ed3dSBarry Smith   PetscFunctionReturn(0);
19439b94acceSBarry Smith }
19449b94acceSBarry Smith 
19454a2ae208SSatish Balay #undef __FUNCT__
19464a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
19479b94acceSBarry Smith /*@C
19483638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19499b94acceSBarry Smith 
1950c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1951c7afd0dbSLois Curfman McInnes 
19529b94acceSBarry Smith    Input Parameter:
19539b94acceSBarry Smith .  snes - the SNES context
19549b94acceSBarry Smith 
19559b94acceSBarry Smith    Output Parameter:
19567bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
195770e92668SMatthew Knepley .  func - the function (or PETSC_NULL)
195870e92668SMatthew Knepley -  ctx - the function context (or PETSC_NULL)
19599b94acceSBarry Smith 
196036851e7fSLois Curfman McInnes    Level: advanced
196136851e7fSLois Curfman McInnes 
1962a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19639b94acceSBarry Smith 
19644b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
19659b94acceSBarry Smith @*/
196670e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
19679b94acceSBarry Smith {
19683a40ed3dSBarry Smith   PetscFunctionBegin;
19694482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19707bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
197100036973SBarry Smith   if (func) *func = snes->computefunction;
197270e92668SMatthew Knepley   if (ctx)  *ctx  = snes->funP;
19733a40ed3dSBarry Smith   PetscFunctionReturn(0);
19749b94acceSBarry Smith }
19759b94acceSBarry Smith 
19764a2ae208SSatish Balay #undef __FUNCT__
19774a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
19783c7409f5SSatish Balay /*@C
19793c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1980d850072dSLois Curfman McInnes    SNES options in the database.
19813c7409f5SSatish Balay 
1982fee21e36SBarry Smith    Collective on SNES
1983fee21e36SBarry Smith 
1984c7afd0dbSLois Curfman McInnes    Input Parameter:
1985c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1986c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
1987c7afd0dbSLois Curfman McInnes 
1988d850072dSLois Curfman McInnes    Notes:
1989a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1990c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
1991d850072dSLois Curfman McInnes 
199236851e7fSLois Curfman McInnes    Level: advanced
199336851e7fSLois Curfman McInnes 
19943c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1995a86d99e1SLois Curfman McInnes 
1996a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
19973c7409f5SSatish Balay @*/
199863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
19993c7409f5SSatish Balay {
2000dfbe8321SBarry Smith   PetscErrorCode ierr;
20013c7409f5SSatish Balay 
20023a40ed3dSBarry Smith   PetscFunctionBegin;
20034482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2004639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
200594b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20063a40ed3dSBarry Smith   PetscFunctionReturn(0);
20073c7409f5SSatish Balay }
20083c7409f5SSatish Balay 
20094a2ae208SSatish Balay #undef __FUNCT__
20104a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
20113c7409f5SSatish Balay /*@C
2012f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2013d850072dSLois Curfman McInnes    SNES options in the database.
20143c7409f5SSatish Balay 
2015fee21e36SBarry Smith    Collective on SNES
2016fee21e36SBarry Smith 
2017c7afd0dbSLois Curfman McInnes    Input Parameters:
2018c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2019c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2020c7afd0dbSLois Curfman McInnes 
2021d850072dSLois Curfman McInnes    Notes:
2022a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2023c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2024d850072dSLois Curfman McInnes 
202536851e7fSLois Curfman McInnes    Level: advanced
202636851e7fSLois Curfman McInnes 
20273c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2028a86d99e1SLois Curfman McInnes 
2029a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20303c7409f5SSatish Balay @*/
203163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
20323c7409f5SSatish Balay {
2033dfbe8321SBarry Smith   PetscErrorCode ierr;
20343c7409f5SSatish Balay 
20353a40ed3dSBarry Smith   PetscFunctionBegin;
20364482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2037639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
203894b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20393a40ed3dSBarry Smith   PetscFunctionReturn(0);
20403c7409f5SSatish Balay }
20413c7409f5SSatish Balay 
20424a2ae208SSatish Balay #undef __FUNCT__
20434a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
20449ab63eb5SSatish Balay /*@C
20453c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20463c7409f5SSatish Balay    SNES options in the database.
20473c7409f5SSatish Balay 
2048c7afd0dbSLois Curfman McInnes    Not Collective
2049c7afd0dbSLois Curfman McInnes 
20503c7409f5SSatish Balay    Input Parameter:
20513c7409f5SSatish Balay .  snes - the SNES context
20523c7409f5SSatish Balay 
20533c7409f5SSatish Balay    Output Parameter:
20543c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20553c7409f5SSatish Balay 
20569ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
20579ab63eb5SSatish Balay    sufficient length to hold the prefix.
20589ab63eb5SSatish Balay 
205936851e7fSLois Curfman McInnes    Level: advanced
206036851e7fSLois Curfman McInnes 
20613c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2062a86d99e1SLois Curfman McInnes 
2063a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20643c7409f5SSatish Balay @*/
2065e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
20663c7409f5SSatish Balay {
2067dfbe8321SBarry Smith   PetscErrorCode ierr;
20683c7409f5SSatish Balay 
20693a40ed3dSBarry Smith   PetscFunctionBegin;
20704482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2071639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
20723a40ed3dSBarry Smith   PetscFunctionReturn(0);
20733c7409f5SSatish Balay }
20743c7409f5SSatish Balay 
2075b2002411SLois Curfman McInnes 
20764a2ae208SSatish Balay #undef __FUNCT__
20774a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
20783cea93caSBarry Smith /*@C
20793cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
20803cea93caSBarry Smith 
20817f6c08e0SMatthew Knepley   Level: advanced
20823cea93caSBarry Smith @*/
208363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2084b2002411SLois Curfman McInnes {
2085e2d1d2b7SBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
2086dfbe8321SBarry Smith   PetscErrorCode ierr;
2087b2002411SLois Curfman McInnes 
2088b2002411SLois Curfman McInnes   PetscFunctionBegin;
2089b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2090c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2091b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2092b2002411SLois Curfman McInnes }
2093da9b6338SBarry Smith 
2094da9b6338SBarry Smith #undef __FUNCT__
2095da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
209663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2097da9b6338SBarry Smith {
2098dfbe8321SBarry Smith   PetscErrorCode ierr;
209977431f27SBarry Smith   PetscInt       N,i,j;
2100da9b6338SBarry Smith   Vec            u,uh,fh;
2101da9b6338SBarry Smith   PetscScalar    value;
2102da9b6338SBarry Smith   PetscReal      norm;
2103da9b6338SBarry Smith 
2104da9b6338SBarry Smith   PetscFunctionBegin;
2105da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2106da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2107da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2108da9b6338SBarry Smith 
2109da9b6338SBarry Smith   /* currently only works for sequential */
2110da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2111da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2112da9b6338SBarry Smith   for (i=0; i<N; i++) {
2113da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
211477431f27SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2115da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2116ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2117da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
21183ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2119da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
212077431f27SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2121da9b6338SBarry Smith       value = -value;
2122da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2123da9b6338SBarry Smith     }
2124da9b6338SBarry Smith   }
2125da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2126da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2127da9b6338SBarry Smith   PetscFunctionReturn(0);
2128da9b6338SBarry Smith }
2129