xref: /petsc/src/snes/interface/snes.c (revision e9a2bbcd8b89561c349212cd9a0aa28814dda9f5)
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
15682738288SBarry Smith .  -snes_monitor - prints residual norm at each iteration
157d132466eSBarry Smith .  -snes_vecmonitor - plots solution at each iteration
158d132466eSBarry Smith .  -snes_vecmonitor_update - plots update to solution at each iteration
15982738288SBarry Smith .  -snes_xmonitor - plots residual norm at each iteration
160e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
1615968eb51SBarry Smith .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
1625968eb51SBarry Smith -  -snes_print_converged_reason - print the reason for convergence/divergence after each solve
16382738288SBarry Smith 
16482738288SBarry Smith     Options Database for Eisenstat-Walker method:
1654b27c08aSLois Curfman McInnes +  -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence
1664b27c08aSLois Curfman McInnes .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
16736851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
16836851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
16936851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
17036851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
17136851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
17236851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
17382738288SBarry Smith 
17411ca99fdSLois Curfman McInnes    Notes:
17511ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
17611ca99fdSLois Curfman McInnes    the users manual.
17783e2fdc7SBarry Smith 
17836851e7fSLois Curfman McInnes    Level: beginner
17936851e7fSLois Curfman McInnes 
1809b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1819b94acceSBarry Smith 
18269ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
1839b94acceSBarry Smith @*/
18463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes)
1859b94acceSBarry Smith {
18694b7f48cSBarry Smith   KSP                 ksp;
187186905e3SBarry Smith   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
188f1af5d2fSBarry Smith   PetscTruth          flg;
189dfbe8321SBarry Smith   PetscErrorCode      ierr;
19077431f27SBarry Smith   PetscInt            i;
1912fc52814SBarry Smith   const char          *deft;
1922fc52814SBarry Smith   char                type[256];
1939b94acceSBarry Smith 
1943a40ed3dSBarry Smith   PetscFunctionBegin;
1954482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
196ca161407SBarry Smith 
197b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
198186905e3SBarry Smith     if (snes->type_name) {
199186905e3SBarry Smith       deft = snes->type_name;
200186905e3SBarry Smith     } else {
2014b27c08aSLois Curfman McInnes       deft = SNESLS;
202d64ed03dSBarry Smith     }
2034bbc92c1SBarry Smith 
204186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
205b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
206d64ed03dSBarry Smith     if (flg) {
207186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
208186905e3SBarry Smith     } else if (!snes->type_name) {
209186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
210d64ed03dSBarry Smith     }
211909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
21293c39befSBarry Smith 
21387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
21470441072SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr);
215186905e3SBarry Smith 
21687828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
217b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
218b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
21950ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
2205968eb51SBarry Smith     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
2215968eb51SBarry Smith     if (flg) {
2225968eb51SBarry Smith       snes->printreason = PETSC_TRUE;
2235968eb51SBarry Smith     }
224186905e3SBarry Smith 
2257aaed0d8SBarry 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);
226186905e3SBarry Smith 
227b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
22887828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
22987828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
23087828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
23187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
23287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
23387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
234186905e3SBarry Smith 
235b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
23693c39befSBarry Smith     if (flg) {snes->converged = 0;}
237b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
2385cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
239b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
240b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
2413a7fca6bSBarry Smith     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
2423a7fca6bSBarry Smith     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
243b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
244b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
245b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
246b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
247b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2487c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
2495ed2d596SBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
2505ed2d596SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
251b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
252186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
253e24b481bSBarry Smith 
254b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2554b27c08aSLois Curfman McInnes     if (flg) {
256186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
25763ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"));CHKERRQ(ierr);
2589b94acceSBarry Smith     }
259639f9d9dSBarry Smith 
26076b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
26176b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
26276b2cf59SMatthew Knepley     }
26376b2cf59SMatthew Knepley 
264186905e3SBarry Smith     if (snes->setfromoptions) {
265186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
266639f9d9dSBarry Smith     }
267639f9d9dSBarry Smith 
268b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2694bbc92c1SBarry Smith 
27094b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
27194b7f48cSBarry Smith   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
27293993e2dSLois Curfman McInnes 
2733a40ed3dSBarry Smith   PetscFunctionReturn(0);
2749b94acceSBarry Smith }
2759b94acceSBarry Smith 
276a847f771SSatish Balay 
2774a2ae208SSatish Balay #undef __FUNCT__
2784a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
2799b94acceSBarry Smith /*@
2809b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2819b94acceSBarry Smith    the nonlinear solvers.
2829b94acceSBarry Smith 
283fee21e36SBarry Smith    Collective on SNES
284fee21e36SBarry Smith 
285c7afd0dbSLois Curfman McInnes    Input Parameters:
286c7afd0dbSLois Curfman McInnes +  snes - the SNES context
287c7afd0dbSLois Curfman McInnes -  usrP - optional user context
288c7afd0dbSLois Curfman McInnes 
28936851e7fSLois Curfman McInnes    Level: intermediate
29036851e7fSLois Curfman McInnes 
2919b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2929b94acceSBarry Smith 
2939b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2949b94acceSBarry Smith @*/
29563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP)
2969b94acceSBarry Smith {
2973a40ed3dSBarry Smith   PetscFunctionBegin;
2984482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2999b94acceSBarry Smith   snes->user		= usrP;
3003a40ed3dSBarry Smith   PetscFunctionReturn(0);
3019b94acceSBarry Smith }
30274679c65SBarry Smith 
3034a2ae208SSatish Balay #undef __FUNCT__
3044a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3059b94acceSBarry Smith /*@C
3069b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3079b94acceSBarry Smith    nonlinear solvers.
3089b94acceSBarry Smith 
309c7afd0dbSLois Curfman McInnes    Not Collective
310c7afd0dbSLois Curfman McInnes 
3119b94acceSBarry Smith    Input Parameter:
3129b94acceSBarry Smith .  snes - SNES context
3139b94acceSBarry Smith 
3149b94acceSBarry Smith    Output Parameter:
3159b94acceSBarry Smith .  usrP - user context
3169b94acceSBarry Smith 
31736851e7fSLois Curfman McInnes    Level: intermediate
31836851e7fSLois Curfman McInnes 
3199b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3209b94acceSBarry Smith 
3219b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3229b94acceSBarry Smith @*/
32363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP)
3249b94acceSBarry Smith {
3253a40ed3dSBarry Smith   PetscFunctionBegin;
3264482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3279b94acceSBarry Smith   *usrP = snes->user;
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
3299b94acceSBarry Smith }
33074679c65SBarry Smith 
3314a2ae208SSatish Balay #undef __FUNCT__
3324a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3339b94acceSBarry Smith /*@
334c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
335c8228a4eSBarry Smith    at this time.
3369b94acceSBarry Smith 
337c7afd0dbSLois Curfman McInnes    Not Collective
338c7afd0dbSLois Curfman McInnes 
3399b94acceSBarry Smith    Input Parameter:
3409b94acceSBarry Smith .  snes - SNES context
3419b94acceSBarry Smith 
3429b94acceSBarry Smith    Output Parameter:
3439b94acceSBarry Smith .  iter - iteration number
3449b94acceSBarry Smith 
345c8228a4eSBarry Smith    Notes:
346c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
347c8228a4eSBarry Smith 
348c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
34908405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
35008405cd6SLois Curfman McInnes .vb
35108405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
35208405cd6SLois Curfman McInnes       if (!(it % 2)) {
35308405cd6SLois Curfman McInnes         [compute Jacobian here]
35408405cd6SLois Curfman McInnes       }
35508405cd6SLois Curfman McInnes .ve
356c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
35708405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
358c8228a4eSBarry Smith 
35936851e7fSLois Curfman McInnes    Level: intermediate
36036851e7fSLois Curfman McInnes 
3619b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3629b94acceSBarry Smith @*/
36363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter)
3649b94acceSBarry Smith {
3653a40ed3dSBarry Smith   PetscFunctionBegin;
3664482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3674482741eSBarry Smith   PetscValidIntPointer(iter,2);
3689b94acceSBarry Smith   *iter = snes->iter;
3693a40ed3dSBarry Smith   PetscFunctionReturn(0);
3709b94acceSBarry Smith }
37174679c65SBarry Smith 
3724a2ae208SSatish Balay #undef __FUNCT__
3734a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
3749b94acceSBarry Smith /*@
3759b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3769b94acceSBarry Smith    with SNESSSetFunction().
3779b94acceSBarry Smith 
378c7afd0dbSLois Curfman McInnes    Collective on SNES
379c7afd0dbSLois Curfman McInnes 
3809b94acceSBarry Smith    Input Parameter:
3819b94acceSBarry Smith .  snes - SNES context
3829b94acceSBarry Smith 
3839b94acceSBarry Smith    Output Parameter:
3849b94acceSBarry Smith .  fnorm - 2-norm of function
3859b94acceSBarry Smith 
38636851e7fSLois Curfman McInnes    Level: intermediate
38736851e7fSLois Curfman McInnes 
3889b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
389a86d99e1SLois Curfman McInnes 
39008405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
3919b94acceSBarry Smith @*/
39263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
3939b94acceSBarry Smith {
3943a40ed3dSBarry Smith   PetscFunctionBegin;
3954482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3964482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
3979b94acceSBarry Smith   *fnorm = snes->norm;
3983a40ed3dSBarry Smith   PetscFunctionReturn(0);
3999b94acceSBarry Smith }
40074679c65SBarry Smith 
4014a2ae208SSatish Balay #undef __FUNCT__
4024a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
4039b94acceSBarry Smith /*@
4049b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4059b94acceSBarry Smith    attempted by the nonlinear solver.
4069b94acceSBarry Smith 
407c7afd0dbSLois Curfman McInnes    Not Collective
408c7afd0dbSLois Curfman McInnes 
4099b94acceSBarry Smith    Input Parameter:
4109b94acceSBarry Smith .  snes - SNES context
4119b94acceSBarry Smith 
4129b94acceSBarry Smith    Output Parameter:
4139b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4149b94acceSBarry Smith 
415c96a6f78SLois Curfman McInnes    Notes:
416c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
417c96a6f78SLois Curfman McInnes 
41836851e7fSLois Curfman McInnes    Level: intermediate
41936851e7fSLois Curfman McInnes 
4209b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4219b94acceSBarry Smith @*/
42263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* nfails)
4239b94acceSBarry Smith {
4243a40ed3dSBarry Smith   PetscFunctionBegin;
4254482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4264482741eSBarry Smith   PetscValidIntPointer(nfails,2);
42750ffb88aSMatthew Knepley   *nfails = snes->numFailures;
42850ffb88aSMatthew Knepley   PetscFunctionReturn(0);
42950ffb88aSMatthew Knepley }
43050ffb88aSMatthew Knepley 
43150ffb88aSMatthew Knepley #undef __FUNCT__
43250ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
43350ffb88aSMatthew Knepley /*@
43450ffb88aSMatthew Knepley    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
43550ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
43650ffb88aSMatthew Knepley 
43750ffb88aSMatthew Knepley    Not Collective
43850ffb88aSMatthew Knepley 
43950ffb88aSMatthew Knepley    Input Parameters:
44050ffb88aSMatthew Knepley +  snes     - SNES context
44150ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
44250ffb88aSMatthew Knepley 
44350ffb88aSMatthew Knepley    Level: intermediate
44450ffb88aSMatthew Knepley 
44550ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
44650ffb88aSMatthew Knepley @*/
44763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt maxFails)
44850ffb88aSMatthew Knepley {
44950ffb88aSMatthew Knepley   PetscFunctionBegin;
4504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
45150ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
45250ffb88aSMatthew Knepley   PetscFunctionReturn(0);
45350ffb88aSMatthew Knepley }
45450ffb88aSMatthew Knepley 
45550ffb88aSMatthew Knepley #undef __FUNCT__
45650ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
45750ffb88aSMatthew Knepley /*@
45850ffb88aSMatthew Knepley    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
45950ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
46050ffb88aSMatthew Knepley 
46150ffb88aSMatthew Knepley    Not Collective
46250ffb88aSMatthew Knepley 
46350ffb88aSMatthew Knepley    Input Parameter:
46450ffb88aSMatthew Knepley .  snes     - SNES context
46550ffb88aSMatthew Knepley 
46650ffb88aSMatthew Knepley    Output Parameter:
46750ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
46850ffb88aSMatthew Knepley 
46950ffb88aSMatthew Knepley    Level: intermediate
47050ffb88aSMatthew Knepley 
47150ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
47250ffb88aSMatthew Knepley @*/
47363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *maxFails)
47450ffb88aSMatthew Knepley {
47550ffb88aSMatthew Knepley   PetscFunctionBegin;
4764482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4774482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
47850ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
4793a40ed3dSBarry Smith   PetscFunctionReturn(0);
4809b94acceSBarry Smith }
481a847f771SSatish Balay 
4824a2ae208SSatish Balay #undef __FUNCT__
4834a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations"
484c96a6f78SLois Curfman McInnes /*@
485c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
486c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
487c96a6f78SLois Curfman McInnes 
488c7afd0dbSLois Curfman McInnes    Not Collective
489c7afd0dbSLois Curfman McInnes 
490c96a6f78SLois Curfman McInnes    Input Parameter:
491c96a6f78SLois Curfman McInnes .  snes - SNES context
492c96a6f78SLois Curfman McInnes 
493c96a6f78SLois Curfman McInnes    Output Parameter:
494c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
495c96a6f78SLois Curfman McInnes 
496c96a6f78SLois Curfman McInnes    Notes:
497c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
498c96a6f78SLois Curfman McInnes 
49936851e7fSLois Curfman McInnes    Level: intermediate
50036851e7fSLois Curfman McInnes 
501c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
502c96a6f78SLois Curfman McInnes @*/
50363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberLinearIterations(SNES snes,PetscInt* lits)
504c96a6f78SLois Curfman McInnes {
5053a40ed3dSBarry Smith   PetscFunctionBegin;
5064482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5074482741eSBarry Smith   PetscValidIntPointer(lits,2);
508c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5093a40ed3dSBarry Smith   PetscFunctionReturn(0);
510c96a6f78SLois Curfman McInnes }
511c96a6f78SLois Curfman McInnes 
5124a2ae208SSatish Balay #undef __FUNCT__
51394b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
5149b94acceSBarry Smith /*@C
51594b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
5169b94acceSBarry Smith 
51794b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
518c7afd0dbSLois Curfman McInnes 
5199b94acceSBarry Smith    Input Parameter:
5209b94acceSBarry Smith .  snes - the SNES context
5219b94acceSBarry Smith 
5229b94acceSBarry Smith    Output Parameter:
52394b7f48cSBarry Smith .  ksp - the KSP context
5249b94acceSBarry Smith 
5259b94acceSBarry Smith    Notes:
52694b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
5279b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5289b94acceSBarry Smith    KSP and PC contexts as well.
5299b94acceSBarry Smith 
53036851e7fSLois Curfman McInnes    Level: beginner
53136851e7fSLois Curfman McInnes 
53294b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
5339b94acceSBarry Smith 
53494b7f48cSBarry Smith .seealso: KSPGetPC()
5359b94acceSBarry Smith @*/
53663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp)
5379b94acceSBarry Smith {
5383a40ed3dSBarry Smith   PetscFunctionBegin;
5394482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5404482741eSBarry Smith   PetscValidPointer(ksp,2);
54194b7f48cSBarry Smith   *ksp = snes->ksp;
5423a40ed3dSBarry Smith   PetscFunctionReturn(0);
5439b94acceSBarry Smith }
54482bf6240SBarry Smith 
5454a2ae208SSatish Balay #undef __FUNCT__
5464a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
5476849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
548e24b481bSBarry Smith {
549e24b481bSBarry Smith   PetscFunctionBegin;
550e24b481bSBarry Smith   PetscFunctionReturn(0);
551e24b481bSBarry Smith }
552e24b481bSBarry Smith 
5539b94acceSBarry Smith /* -----------------------------------------------------------*/
5544a2ae208SSatish Balay #undef __FUNCT__
5554a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
5569b94acceSBarry Smith /*@C
5579b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5589b94acceSBarry Smith 
559c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
560c7afd0dbSLois Curfman McInnes 
561c7afd0dbSLois Curfman McInnes    Input Parameters:
562c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
5639b94acceSBarry Smith 
5649b94acceSBarry Smith    Output Parameter:
5659b94acceSBarry Smith .  outsnes - the new SNES context
5669b94acceSBarry Smith 
567c7afd0dbSLois Curfman McInnes    Options Database Keys:
568c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
569c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
570c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
571c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
572c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
573c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
574c1f60f51SBarry Smith 
57536851e7fSLois Curfman McInnes    Level: beginner
57636851e7fSLois Curfman McInnes 
5779b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5789b94acceSBarry Smith 
5794b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES
5809b94acceSBarry Smith @*/
58163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes)
5829b94acceSBarry Smith {
583dfbe8321SBarry Smith   PetscErrorCode      ierr;
5849b94acceSBarry Smith   SNES                snes;
5859b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
58637fcc0dbSBarry Smith 
5873a40ed3dSBarry Smith   PetscFunctionBegin;
588ed1caa07SMatthew Knepley   PetscValidPointer(outsnes,2);
5898ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
5908ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
5918ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
5928ba1e511SMatthew Knepley #endif
5938ba1e511SMatthew Knepley 
59452e6d16bSBarry Smith   ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr);
595e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
5969b94acceSBarry Smith   snes->max_its           = 50;
5979750a799SBarry Smith   snes->max_funcs	  = 10000;
5989b94acceSBarry Smith   snes->norm		  = 0.0;
599b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
600b4874afaSBarry Smith   snes->ttol              = 0.0;
60170441072SBarry Smith   snes->abstol		  = 1.e-50;
6029b94acceSBarry Smith   snes->xtol		  = 1.e-8;
6034b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
6049b94acceSBarry Smith   snes->nfuncs            = 0;
60550ffb88aSMatthew Knepley   snes->numFailures       = 0;
60650ffb88aSMatthew Knepley   snes->maxFailures       = 1;
6077a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
608639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6099b94acceSBarry Smith   snes->data              = 0;
6109b94acceSBarry Smith   snes->view              = 0;
6114dc4c822SBarry Smith   snes->setupcalled       = PETSC_FALSE;
612186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
6136f24a144SLois Curfman McInnes   snes->vwork             = 0;
6146f24a144SLois Curfman McInnes   snes->nwork             = 0;
615758f92a0SBarry Smith   snes->conv_hist_len     = 0;
616758f92a0SBarry Smith   snes->conv_hist_max     = 0;
617758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
618758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
619758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
620184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
6219b94acceSBarry Smith 
6229b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
623b0a32e0cSBarry Smith   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
62452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));CHKERRQ(ierr);
6259b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6269b94acceSBarry Smith   kctx->version     = 2;
6279b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6289b94acceSBarry Smith                              this was too large for some test cases */
6299b94acceSBarry Smith   kctx->rtol_last   = 0;
6309b94acceSBarry Smith   kctx->rtol_max    = .9;
6319b94acceSBarry Smith   kctx->gamma       = 1.0;
6329b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6339b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6349b94acceSBarry Smith   kctx->threshold   = .1;
6359b94acceSBarry Smith   kctx->lresid_last = 0;
6369b94acceSBarry Smith   kctx->norm_last   = 0;
6379b94acceSBarry Smith 
63894b7f48cSBarry Smith   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
63952e6d16bSBarry Smith   ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr);
6405334005bSBarry Smith 
6419b94acceSBarry Smith   *outsnes = snes;
64200036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
6433a40ed3dSBarry Smith   PetscFunctionReturn(0);
6449b94acceSBarry Smith }
6459b94acceSBarry Smith 
6464a2ae208SSatish Balay #undef __FUNCT__
6474a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
6489b94acceSBarry Smith /*@C
6499b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6509b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6519b94acceSBarry Smith    equations.
6529b94acceSBarry Smith 
653fee21e36SBarry Smith    Collective on SNES
654fee21e36SBarry Smith 
655c7afd0dbSLois Curfman McInnes    Input Parameters:
656c7afd0dbSLois Curfman McInnes +  snes - the SNES context
657c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
658c7afd0dbSLois Curfman McInnes .  r - vector to store function value
659c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
660c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6619b94acceSBarry Smith 
662c7afd0dbSLois Curfman McInnes    Calling sequence of func:
6638d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
664c7afd0dbSLois Curfman McInnes 
665313e4042SLois Curfman McInnes .  f - function vector
666c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6679b94acceSBarry Smith 
6689b94acceSBarry Smith    Notes:
6699b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6709b94acceSBarry Smith $      f'(x) x = -f(x),
671c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6729b94acceSBarry Smith 
67336851e7fSLois Curfman McInnes    Level: beginner
67436851e7fSLois Curfman McInnes 
6759b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6769b94acceSBarry Smith 
677a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6789b94acceSBarry Smith @*/
67963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
6809b94acceSBarry Smith {
6813a40ed3dSBarry Smith   PetscFunctionBegin;
6824482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
6834482741eSBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
684c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,r,2);
685184914b5SBarry Smith 
6869b94acceSBarry Smith   snes->computefunction     = func;
6879b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6889b94acceSBarry Smith   snes->funP                = ctx;
6893a40ed3dSBarry Smith   PetscFunctionReturn(0);
6909b94acceSBarry Smith }
6919b94acceSBarry Smith 
6923ab0aad5SBarry Smith /* --------------------------------------------------------------- */
6933ab0aad5SBarry Smith #undef __FUNCT__
6943ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs"
6953ab0aad5SBarry Smith /*@C
6963ab0aad5SBarry Smith    SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set
6973ab0aad5SBarry Smith    it assumes a zero right hand side.
6983ab0aad5SBarry Smith 
6993ab0aad5SBarry Smith    Collective on SNES
7003ab0aad5SBarry Smith 
7013ab0aad5SBarry Smith    Input Parameters:
7023ab0aad5SBarry Smith +  snes - the SNES context
7033ab0aad5SBarry Smith -  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7043ab0aad5SBarry Smith 
7053ab0aad5SBarry Smith    Level: intermediate
7063ab0aad5SBarry Smith 
7073ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side
7083ab0aad5SBarry Smith 
7091096aae1SMatthew Knepley .seealso: SNESGetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7103ab0aad5SBarry Smith @*/
71163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetRhs(SNES snes,Vec rhs)
7123ab0aad5SBarry Smith {
713dfbe8321SBarry Smith   PetscErrorCode ierr;
7143ab0aad5SBarry Smith 
7153ab0aad5SBarry Smith   PetscFunctionBegin;
7163ab0aad5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7173ab0aad5SBarry Smith   if (rhs) {
7183ab0aad5SBarry Smith     PetscValidHeaderSpecific(rhs,VEC_COOKIE,2);
7193ab0aad5SBarry Smith     PetscCheckSameComm(snes,1,rhs,2);
7203ab0aad5SBarry Smith     ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr);
7213ab0aad5SBarry Smith   }
7223ab0aad5SBarry Smith   if (snes->afine) {
7233ab0aad5SBarry Smith     ierr = VecDestroy(snes->afine);CHKERRQ(ierr);
7243ab0aad5SBarry Smith   }
7253ab0aad5SBarry Smith   snes->afine = rhs;
7263ab0aad5SBarry Smith   PetscFunctionReturn(0);
7273ab0aad5SBarry Smith }
7283ab0aad5SBarry Smith 
7294a2ae208SSatish Balay #undef __FUNCT__
7301096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs"
7311096aae1SMatthew Knepley /*@C
7321096aae1SMatthew Knepley    SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set
7331096aae1SMatthew Knepley    it assumes a zero right hand side.
7341096aae1SMatthew Knepley 
7351096aae1SMatthew Knepley    Collective on SNES
7361096aae1SMatthew Knepley 
7371096aae1SMatthew Knepley    Input Parameter:
7381096aae1SMatthew Knepley .  snes - the SNES context
7391096aae1SMatthew Knepley 
7401096aae1SMatthew Knepley    Output Parameter:
7411096aae1SMatthew Knepley .  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7421096aae1SMatthew Knepley 
7431096aae1SMatthew Knepley    Level: intermediate
7441096aae1SMatthew Knepley 
7451096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side
7461096aae1SMatthew Knepley 
7471096aae1SMatthew Knepley .seealso: SNESSetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7481096aae1SMatthew Knepley @*/
7491096aae1SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs)
7501096aae1SMatthew Knepley {
7511096aae1SMatthew Knepley   PetscFunctionBegin;
7521096aae1SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7531096aae1SMatthew Knepley   PetscValidPointer(rhs,2);
7541096aae1SMatthew Knepley   *rhs = snes->afine;
7551096aae1SMatthew Knepley   PetscFunctionReturn(0);
7561096aae1SMatthew Knepley }
7571096aae1SMatthew Knepley 
7581096aae1SMatthew Knepley #undef __FUNCT__
7594a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
7609b94acceSBarry Smith /*@
76136851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
7629b94acceSBarry Smith                          SNESSetFunction().
7639b94acceSBarry Smith 
764c7afd0dbSLois Curfman McInnes    Collective on SNES
765c7afd0dbSLois Curfman McInnes 
7669b94acceSBarry Smith    Input Parameters:
767c7afd0dbSLois Curfman McInnes +  snes - the SNES context
768c7afd0dbSLois Curfman McInnes -  x - input vector
7699b94acceSBarry Smith 
7709b94acceSBarry Smith    Output Parameter:
7713638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7729b94acceSBarry Smith 
7731bffabb2SLois Curfman McInnes    Notes:
77436851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
77536851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
77636851e7fSLois Curfman McInnes    themselves.
77736851e7fSLois Curfman McInnes 
77836851e7fSLois Curfman McInnes    Level: developer
77936851e7fSLois Curfman McInnes 
7809b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7819b94acceSBarry Smith 
782a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7839b94acceSBarry Smith @*/
78463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y)
7859b94acceSBarry Smith {
786dfbe8321SBarry Smith   PetscErrorCode ierr;
7879b94acceSBarry Smith 
7883a40ed3dSBarry Smith   PetscFunctionBegin;
7894482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7904482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
7914482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
792c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
793c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
794184914b5SBarry Smith 
795d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
7961096aae1SMatthew Knepley   if (snes->computefunction) {
797d64ed03dSBarry Smith     PetscStackPush("SNES user function");
798*e9a2bbcdSBarry Smith     CHKMEMQ;
79919717074SBarry Smith     ierr = (*snes->computefunction)(snes,x,y,snes->funP);
800*e9a2bbcdSBarry Smith     CHKMEMQ;
801d64ed03dSBarry Smith     PetscStackPop;
802d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
80319717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr);
80419717074SBarry Smith     }
805d5e45103SBarry Smith     CHKERRQ(ierr);
8061096aae1SMatthew Knepley   } else if (snes->afine) {
8071096aae1SMatthew Knepley     ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr);
8081096aae1SMatthew Knepley   } else {
8091096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve().");
8101096aae1SMatthew Knepley   }
8113ab0aad5SBarry Smith   if (snes->afine) {
8123ab0aad5SBarry Smith     PetscScalar mone = -1.0;
8132dcb1b2aSMatthew Knepley     ierr = VecAXPY(y,mone,snes->afine);CHKERRQ(ierr);
8143ab0aad5SBarry Smith   }
815ae3c334cSLois Curfman McInnes   snes->nfuncs++;
816d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8173a40ed3dSBarry Smith   PetscFunctionReturn(0);
8189b94acceSBarry Smith }
8199b94acceSBarry Smith 
8204a2ae208SSatish Balay #undef __FUNCT__
8214a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
82262fef451SLois Curfman McInnes /*@
82362fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
82462fef451SLois Curfman McInnes    set with SNESSetJacobian().
82562fef451SLois Curfman McInnes 
826c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
827c7afd0dbSLois Curfman McInnes 
82862fef451SLois Curfman McInnes    Input Parameters:
829c7afd0dbSLois Curfman McInnes +  snes - the SNES context
830c7afd0dbSLois Curfman McInnes -  x - input vector
83162fef451SLois Curfman McInnes 
83262fef451SLois Curfman McInnes    Output Parameters:
833c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
83462fef451SLois Curfman McInnes .  B - optional preconditioning matrix
835c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
836fee21e36SBarry Smith 
83762fef451SLois Curfman McInnes    Notes:
83862fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83962fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
84062fef451SLois Curfman McInnes 
84194b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
842dc5a77f8SLois Curfman McInnes    flag parameter.
84362fef451SLois Curfman McInnes 
84436851e7fSLois Curfman McInnes    Level: developer
84536851e7fSLois Curfman McInnes 
84662fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
84762fef451SLois Curfman McInnes 
84894b7f48cSBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators()
84962fef451SLois Curfman McInnes @*/
85063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8519b94acceSBarry Smith {
852dfbe8321SBarry Smith   PetscErrorCode ierr;
8533a40ed3dSBarry Smith 
8543a40ed3dSBarry Smith   PetscFunctionBegin;
8554482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8564482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
8574482741eSBarry Smith   PetscValidPointer(flg,5);
858c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
8593a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
860d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
861c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
862d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
8639b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
864d64ed03dSBarry Smith   PetscStackPop;
865d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
8666d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
8674482741eSBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
8684482741eSBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE,4);
8693a40ed3dSBarry Smith   PetscFunctionReturn(0);
8709b94acceSBarry Smith }
8719b94acceSBarry Smith 
8724a2ae208SSatish Balay #undef __FUNCT__
8734a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
8749b94acceSBarry Smith /*@C
8759b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
876044dda88SLois Curfman McInnes    location to store the matrix.
8779b94acceSBarry Smith 
878c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
879c7afd0dbSLois Curfman McInnes 
8809b94acceSBarry Smith    Input Parameters:
881c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8829b94acceSBarry Smith .  A - Jacobian matrix
8839b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8849b94acceSBarry Smith .  func - Jacobian evaluation routine
885c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
8862cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
8879b94acceSBarry Smith 
8889b94acceSBarry Smith    Calling sequence of func:
8898d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
8909b94acceSBarry Smith 
891c7afd0dbSLois Curfman McInnes +  x - input vector
8929b94acceSBarry Smith .  A - Jacobian matrix
8939b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
894ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
89594b7f48cSBarry Smith    structure (same as flag in KSPSetOperators())
896c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
8979b94acceSBarry Smith 
8989b94acceSBarry Smith    Notes:
89994b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
9002cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
901ac21db08SLois Curfman McInnes 
902ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9039b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9049b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9059b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9069b94acceSBarry Smith    throughout the global iterations.
9079b94acceSBarry Smith 
90836851e7fSLois Curfman McInnes    Level: beginner
90936851e7fSLois Curfman McInnes 
9109b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9119b94acceSBarry Smith 
91213f74950SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), , MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor()
9139b94acceSBarry Smith @*/
91463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
9159b94acceSBarry Smith {
916dfbe8321SBarry Smith   PetscErrorCode ierr;
9173a7fca6bSBarry Smith 
9183a40ed3dSBarry Smith   PetscFunctionBegin;
9194482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9204482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
9214482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
922c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
923c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
9243a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
9253a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
9263a7fca6bSBarry Smith   if (A) {
9273a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
9289b94acceSBarry Smith     snes->jacobian = A;
9293a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
9303a7fca6bSBarry Smith   }
9313a7fca6bSBarry Smith   if (B) {
9323a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
9339b94acceSBarry Smith     snes->jacobian_pre = B;
9343a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
9353a7fca6bSBarry Smith   }
93669a612a9SBarry Smith   ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9373a40ed3dSBarry Smith   PetscFunctionReturn(0);
9389b94acceSBarry Smith }
93962fef451SLois Curfman McInnes 
9404a2ae208SSatish Balay #undef __FUNCT__
9414a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
942c2aafc4cSSatish Balay /*@C
943b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
944b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
945b4fd4287SBarry Smith 
946c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
947c7afd0dbSLois Curfman McInnes 
948b4fd4287SBarry Smith    Input Parameter:
949b4fd4287SBarry Smith .  snes - the nonlinear solver context
950b4fd4287SBarry Smith 
951b4fd4287SBarry Smith    Output Parameters:
952c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
953b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
95470e92668SMatthew Knepley .  func - location to put Jacobian function (or PETSC_NULL)
95570e92668SMatthew Knepley -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
956fee21e36SBarry Smith 
95736851e7fSLois Curfman McInnes    Level: advanced
95836851e7fSLois Curfman McInnes 
959b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
960b4fd4287SBarry Smith @*/
96170e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
962b4fd4287SBarry Smith {
9633a40ed3dSBarry Smith   PetscFunctionBegin;
9644482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
965b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
966b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
96700036973SBarry Smith   if (func) *func = snes->computejacobian;
96870e92668SMatthew Knepley   if (ctx)  *ctx  = snes->jacP;
9693a40ed3dSBarry Smith   PetscFunctionReturn(0);
970b4fd4287SBarry Smith }
971b4fd4287SBarry Smith 
9729b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
97363dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
9749b94acceSBarry Smith 
9754a2ae208SSatish Balay #undef __FUNCT__
9764a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
9779b94acceSBarry Smith /*@
9789b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
979272ac6f2SLois Curfman McInnes    of a nonlinear solver.
9809b94acceSBarry Smith 
981fee21e36SBarry Smith    Collective on SNES
982fee21e36SBarry Smith 
983c7afd0dbSLois Curfman McInnes    Input Parameters:
98470e92668SMatthew Knepley .  snes - the SNES context
985c7afd0dbSLois Curfman McInnes 
986272ac6f2SLois Curfman McInnes    Notes:
987272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
988272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
989272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
990272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
991272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
992272ac6f2SLois Curfman McInnes 
99336851e7fSLois Curfman McInnes    Level: advanced
99436851e7fSLois Curfman McInnes 
9959b94acceSBarry Smith .keywords: SNES, nonlinear, setup
9969b94acceSBarry Smith 
9979b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
9989b94acceSBarry Smith @*/
99970e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
10009b94acceSBarry Smith {
1001dfbe8321SBarry Smith   PetscErrorCode ierr;
10024b27c08aSLois Curfman McInnes   PetscTruth     flg, iseqtr;
10033a40ed3dSBarry Smith 
10043a40ed3dSBarry Smith   PetscFunctionBegin;
10054482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10064dc4c822SBarry Smith   if (snes->setupcalled) PetscFunctionReturn(0);
10079b94acceSBarry Smith 
1008b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1009c1f60f51SBarry Smith   /*
1010c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1011dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1012c1f60f51SBarry Smith   */
1013c1f60f51SBarry Smith   if (flg) {
1014c1f60f51SBarry Smith     Mat J;
10155a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10165a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
101763ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr);
10183a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
10193a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1020c1f60f51SBarry Smith   }
102145fc7adcSBarry Smith 
1022c8d321e0SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE)
102345fc7adcSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
102445fc7adcSBarry Smith   if (flg) {
102545fc7adcSBarry Smith     Mat J;
102645fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
102745fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
102845fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
102945fc7adcSBarry Smith   }
103032a4b47aSMatthew Knepley #endif
103145fc7adcSBarry Smith 
1032b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1033c1f60f51SBarry Smith   /*
1034dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1035c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1036c1f60f51SBarry Smith    */
103731615d04SLois Curfman McInnes   if (flg) {
1038272ac6f2SLois Curfman McInnes     Mat  J;
1039b5d62d44SBarry Smith     KSP ksp;
104094b7f48cSBarry Smith     PC   pc;
1041f3ef73ceSBarry Smith 
10425a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10433a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
104463ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr);
10453a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
10463a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
10473a7fca6bSBarry Smith 
1048f3ef73ceSBarry Smith     /* force no preconditioner */
104994b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1050b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1051a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1052a9815358SBarry Smith     if (!flg) {
1053f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1054272ac6f2SLois Curfman McInnes     }
1055a9815358SBarry Smith   }
1056f3ef73ceSBarry Smith 
10571096aae1SMatthew Knepley   if (!snes->vec_func && !snes->afine) {
10581096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10591096aae1SMatthew Knepley   }
10601096aae1SMatthew Knepley   if (!snes->computefunction && !snes->afine) {
10611096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10621096aae1SMatthew Knepley   }
106329bbc08cSBarry Smith   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1064a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
106529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1066a8c6a408SBarry Smith   }
1067a703fe33SLois Curfman McInnes 
1068a703fe33SLois Curfman McInnes   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
10694b27c08aSLois Curfman McInnes   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
10706831982aSBarry Smith   if (snes->ksp_ewconv && !iseqtr) {
107194b7f48cSBarry Smith     KSP ksp;
107294b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1073186905e3SBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1074a703fe33SLois Curfman McInnes   }
10754b27c08aSLois Curfman McInnes 
1076a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
10777aaed0d8SBarry Smith   snes->setupcalled = PETSC_TRUE;
10783a40ed3dSBarry Smith   PetscFunctionReturn(0);
10799b94acceSBarry Smith }
10809b94acceSBarry Smith 
10814a2ae208SSatish Balay #undef __FUNCT__
10824a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
10839b94acceSBarry Smith /*@C
10849b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
10859b94acceSBarry Smith    with SNESCreate().
10869b94acceSBarry Smith 
1087c7afd0dbSLois Curfman McInnes    Collective on SNES
1088c7afd0dbSLois Curfman McInnes 
10899b94acceSBarry Smith    Input Parameter:
10909b94acceSBarry Smith .  snes - the SNES context
10919b94acceSBarry Smith 
109236851e7fSLois Curfman McInnes    Level: beginner
109336851e7fSLois Curfman McInnes 
10949b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
10959b94acceSBarry Smith 
109663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
10979b94acceSBarry Smith @*/
109863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
10999b94acceSBarry Smith {
110077431f27SBarry Smith   PetscInt       i;
11016849ba73SBarry Smith   PetscErrorCode ierr;
11023a40ed3dSBarry Smith 
11033a40ed3dSBarry Smith   PetscFunctionBegin;
11044482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
11053a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1106d4bb536fSBarry Smith 
1107be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
11080f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1109be0abb6dSBarry Smith 
1110e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1111606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
11123a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
11133a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
11143ab0aad5SBarry Smith   if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);}
111594b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1116522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1117b8a78c4aSBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1118b8a78c4aSBarry Smith     if (snes->monitordestroy[i]) {
1119b8a78c4aSBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1120b8a78c4aSBarry Smith     }
1121b8a78c4aSBarry Smith   }
1122a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
11233a40ed3dSBarry Smith   PetscFunctionReturn(0);
11249b94acceSBarry Smith }
11259b94acceSBarry Smith 
11269b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11279b94acceSBarry Smith 
11284a2ae208SSatish Balay #undef __FUNCT__
11294a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
11309b94acceSBarry Smith /*@
1131d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11329b94acceSBarry Smith 
1133c7afd0dbSLois Curfman McInnes    Collective on SNES
1134c7afd0dbSLois Curfman McInnes 
11359b94acceSBarry Smith    Input Parameters:
1136c7afd0dbSLois Curfman McInnes +  snes - the SNES context
113770441072SBarry Smith .  abstol - absolute convergence tolerance
113833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
113933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
114033174efeSLois Curfman McInnes            of the change in the solution between steps
114133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1142c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1143fee21e36SBarry Smith 
114433174efeSLois Curfman McInnes    Options Database Keys:
114570441072SBarry Smith +    -snes_atol <abstol> - Sets abstol
1146c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1147c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1148c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1149c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
11509b94acceSBarry Smith 
1151d7a720efSLois Curfman McInnes    Notes:
11529b94acceSBarry Smith    The default maximum number of iterations is 50.
11539b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11549b94acceSBarry Smith 
115536851e7fSLois Curfman McInnes    Level: intermediate
115636851e7fSLois Curfman McInnes 
115733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11589b94acceSBarry Smith 
11592492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance()
11609b94acceSBarry Smith @*/
116163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
11629b94acceSBarry Smith {
11633a40ed3dSBarry Smith   PetscFunctionBegin;
11644482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
116570441072SBarry Smith   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1166d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1167d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1168d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1169d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11703a40ed3dSBarry Smith   PetscFunctionReturn(0);
11719b94acceSBarry Smith }
11729b94acceSBarry Smith 
11734a2ae208SSatish Balay #undef __FUNCT__
11744a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
11759b94acceSBarry Smith /*@
117633174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
117733174efeSLois Curfman McInnes 
1178c7afd0dbSLois Curfman McInnes    Not Collective
1179c7afd0dbSLois Curfman McInnes 
118033174efeSLois Curfman McInnes    Input Parameters:
1181c7afd0dbSLois Curfman McInnes +  snes - the SNES context
118270441072SBarry Smith .  abstol - absolute convergence tolerance
118333174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
118433174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
118533174efeSLois Curfman McInnes            of the change in the solution between steps
118633174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1187c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1188fee21e36SBarry Smith 
118933174efeSLois Curfman McInnes    Notes:
119033174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
119133174efeSLois Curfman McInnes 
119236851e7fSLois Curfman McInnes    Level: intermediate
119336851e7fSLois Curfman McInnes 
119433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
119533174efeSLois Curfman McInnes 
119633174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
119733174efeSLois Curfman McInnes @*/
119863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
119933174efeSLois Curfman McInnes {
12003a40ed3dSBarry Smith   PetscFunctionBegin;
12014482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
120270441072SBarry Smith   if (abstol)  *abstol  = snes->abstol;
120333174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
120433174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
120533174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
120633174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
12073a40ed3dSBarry Smith   PetscFunctionReturn(0);
120833174efeSLois Curfman McInnes }
120933174efeSLois Curfman McInnes 
12104a2ae208SSatish Balay #undef __FUNCT__
12114a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
121233174efeSLois Curfman McInnes /*@
12139b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12149b94acceSBarry Smith 
1215fee21e36SBarry Smith    Collective on SNES
1216fee21e36SBarry Smith 
1217c7afd0dbSLois Curfman McInnes    Input Parameters:
1218c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1219c7afd0dbSLois Curfman McInnes -  tol - tolerance
1220c7afd0dbSLois Curfman McInnes 
12219b94acceSBarry Smith    Options Database Key:
1222c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
12239b94acceSBarry Smith 
122436851e7fSLois Curfman McInnes    Level: intermediate
122536851e7fSLois Curfman McInnes 
12269b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12279b94acceSBarry Smith 
12282492ecdbSBarry Smith .seealso: SNESSetTolerances()
12299b94acceSBarry Smith @*/
123063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
12319b94acceSBarry Smith {
12323a40ed3dSBarry Smith   PetscFunctionBegin;
12334482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12349b94acceSBarry Smith   snes->deltatol = tol;
12353a40ed3dSBarry Smith   PetscFunctionReturn(0);
12369b94acceSBarry Smith }
12379b94acceSBarry Smith 
1238df9fa365SBarry Smith /*
1239df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1240df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1241df9fa365SBarry Smith    macros instead of functions
1242df9fa365SBarry Smith */
12434a2ae208SSatish Balay #undef __FUNCT__
12444a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
124563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1246ce1608b8SBarry Smith {
1247dfbe8321SBarry Smith   PetscErrorCode ierr;
1248ce1608b8SBarry Smith 
1249ce1608b8SBarry Smith   PetscFunctionBegin;
12504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1251ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1252ce1608b8SBarry Smith   PetscFunctionReturn(0);
1253ce1608b8SBarry Smith }
1254ce1608b8SBarry Smith 
12554a2ae208SSatish Balay #undef __FUNCT__
12564a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
125763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1258df9fa365SBarry Smith {
1259dfbe8321SBarry Smith   PetscErrorCode ierr;
1260df9fa365SBarry Smith 
1261df9fa365SBarry Smith   PetscFunctionBegin;
1262df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1263df9fa365SBarry Smith   PetscFunctionReturn(0);
1264df9fa365SBarry Smith }
1265df9fa365SBarry Smith 
12664a2ae208SSatish Balay #undef __FUNCT__
12674a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
126863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw)
1269df9fa365SBarry Smith {
1270dfbe8321SBarry Smith   PetscErrorCode ierr;
1271df9fa365SBarry Smith 
1272df9fa365SBarry Smith   PetscFunctionBegin;
1273df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1274df9fa365SBarry Smith   PetscFunctionReturn(0);
1275df9fa365SBarry Smith }
1276df9fa365SBarry Smith 
12779b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12789b94acceSBarry Smith 
12794a2ae208SSatish Balay #undef __FUNCT__
12804a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
12819b94acceSBarry Smith /*@C
12825cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
12839b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12849b94acceSBarry Smith    progress.
12859b94acceSBarry Smith 
1286fee21e36SBarry Smith    Collective on SNES
1287fee21e36SBarry Smith 
1288c7afd0dbSLois Curfman McInnes    Input Parameters:
1289c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1290c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1291b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1292b3006f0bSLois Curfman McInnes           monitor routine (use PETSC_NULL if no context is desitre)
1293b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1294b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
12959b94acceSBarry Smith 
1296c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1297a7cc72afSBarry Smith $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1298c7afd0dbSLois Curfman McInnes 
1299c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1300c7afd0dbSLois Curfman McInnes .    its - iteration number
1301c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
130240a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
13039b94acceSBarry Smith 
13049665c990SLois Curfman McInnes    Options Database Keys:
1305c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1306c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1307c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1308c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1309c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1310c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1311c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1312c7afd0dbSLois Curfman McInnes                             the options database.
13139665c990SLois Curfman McInnes 
1314639f9d9dSBarry Smith    Notes:
13156bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13166bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13176bc08f3fSLois Curfman McInnes    order in which they were set.
1318639f9d9dSBarry Smith 
131936851e7fSLois Curfman McInnes    Level: intermediate
132036851e7fSLois Curfman McInnes 
13219b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13229b94acceSBarry Smith 
13235cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
13249b94acceSBarry Smith @*/
132563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
13269b94acceSBarry Smith {
13273a40ed3dSBarry Smith   PetscFunctionBegin;
13284482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1329639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
133029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1331639f9d9dSBarry Smith   }
1332639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1333b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1334639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13353a40ed3dSBarry Smith   PetscFunctionReturn(0);
13369b94acceSBarry Smith }
13379b94acceSBarry Smith 
13384a2ae208SSatish Balay #undef __FUNCT__
13394a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
13405cd90555SBarry Smith /*@C
13415cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
13425cd90555SBarry Smith 
1343c7afd0dbSLois Curfman McInnes    Collective on SNES
1344c7afd0dbSLois Curfman McInnes 
13455cd90555SBarry Smith    Input Parameters:
13465cd90555SBarry Smith .  snes - the SNES context
13475cd90555SBarry Smith 
13481a480d89SAdministrator    Options Database Key:
1349c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1350c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1351c7afd0dbSLois Curfman McInnes     set via the options database
13525cd90555SBarry Smith 
13535cd90555SBarry Smith    Notes:
13545cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
13555cd90555SBarry Smith 
135636851e7fSLois Curfman McInnes    Level: intermediate
135736851e7fSLois Curfman McInnes 
13585cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
13595cd90555SBarry Smith 
13605cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
13615cd90555SBarry Smith @*/
136263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes)
13635cd90555SBarry Smith {
13645cd90555SBarry Smith   PetscFunctionBegin;
13654482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13665cd90555SBarry Smith   snes->numbermonitors = 0;
13675cd90555SBarry Smith   PetscFunctionReturn(0);
13685cd90555SBarry Smith }
13695cd90555SBarry Smith 
13704a2ae208SSatish Balay #undef __FUNCT__
13714a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
13729b94acceSBarry Smith /*@C
13739b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13749b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13759b94acceSBarry Smith 
1376fee21e36SBarry Smith    Collective on SNES
1377fee21e36SBarry Smith 
1378c7afd0dbSLois Curfman McInnes    Input Parameters:
1379c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1380c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1381c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1382c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
13839b94acceSBarry Smith 
1384c7afd0dbSLois Curfman McInnes    Calling sequence of func:
13856849ba73SBarry Smith $     PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1386c7afd0dbSLois Curfman McInnes 
1387c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1388c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1389184914b5SBarry Smith .    reason - reason for convergence/divergence
1390c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
13914b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
13924b27c08aSLois Curfman McInnes -    f - 2-norm of function
13939b94acceSBarry Smith 
139436851e7fSLois Curfman McInnes    Level: advanced
139536851e7fSLois Curfman McInnes 
13969b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
13979b94acceSBarry Smith 
13984b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR()
13999b94acceSBarry Smith @*/
140063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
14019b94acceSBarry Smith {
14023a40ed3dSBarry Smith   PetscFunctionBegin;
14034482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14049b94acceSBarry Smith   (snes)->converged = func;
14059b94acceSBarry Smith   (snes)->cnvP      = cctx;
14063a40ed3dSBarry Smith   PetscFunctionReturn(0);
14079b94acceSBarry Smith }
14089b94acceSBarry Smith 
14094a2ae208SSatish Balay #undef __FUNCT__
14104a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
1411184914b5SBarry Smith /*@C
1412184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1413184914b5SBarry Smith 
1414184914b5SBarry Smith    Not Collective
1415184914b5SBarry Smith 
1416184914b5SBarry Smith    Input Parameter:
1417184914b5SBarry Smith .  snes - the SNES context
1418184914b5SBarry Smith 
1419184914b5SBarry Smith    Output Parameter:
1420e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1421184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1422184914b5SBarry Smith 
1423184914b5SBarry Smith    Level: intermediate
1424184914b5SBarry Smith 
1425184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1426184914b5SBarry Smith 
1427184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1428184914b5SBarry Smith 
14294b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1430184914b5SBarry Smith @*/
143163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1432184914b5SBarry Smith {
1433184914b5SBarry Smith   PetscFunctionBegin;
14344482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14354482741eSBarry Smith   PetscValidPointer(reason,2);
1436184914b5SBarry Smith   *reason = snes->reason;
1437184914b5SBarry Smith   PetscFunctionReturn(0);
1438184914b5SBarry Smith }
1439184914b5SBarry Smith 
14404a2ae208SSatish Balay #undef __FUNCT__
14414a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1442c9005455SLois Curfman McInnes /*@
1443c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1444c9005455SLois Curfman McInnes 
1445fee21e36SBarry Smith    Collective on SNES
1446fee21e36SBarry Smith 
1447c7afd0dbSLois Curfman McInnes    Input Parameters:
1448c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1449c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1450cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1451758f92a0SBarry Smith .  na  - size of a and its
145264731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1453758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1454c7afd0dbSLois Curfman McInnes 
1455c9005455SLois Curfman McInnes    Notes:
14564b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1457c9005455SLois Curfman McInnes    at each step.
1458c9005455SLois Curfman McInnes 
1459c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1460c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1461c9005455SLois Curfman McInnes    during the section of code that is being timed.
1462c9005455SLois Curfman McInnes 
146336851e7fSLois Curfman McInnes    Level: intermediate
146436851e7fSLois Curfman McInnes 
1465c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1466758f92a0SBarry Smith 
146708405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1468758f92a0SBarry Smith 
1469c9005455SLois Curfman McInnes @*/
147063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset)
1471c9005455SLois Curfman McInnes {
14723a40ed3dSBarry Smith   PetscFunctionBegin;
14734482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14744482741eSBarry Smith   if (na) PetscValidScalarPointer(a,2);
1475c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1476758f92a0SBarry Smith   snes->conv_hist_its   = its;
1477758f92a0SBarry Smith   snes->conv_hist_max   = na;
1478758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1479758f92a0SBarry Smith   PetscFunctionReturn(0);
1480758f92a0SBarry Smith }
1481758f92a0SBarry Smith 
14824a2ae208SSatish Balay #undef __FUNCT__
14834a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
14840c4c9dddSBarry Smith /*@C
1485758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1486758f92a0SBarry Smith 
1487758f92a0SBarry Smith    Collective on SNES
1488758f92a0SBarry Smith 
1489758f92a0SBarry Smith    Input Parameter:
1490758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1491758f92a0SBarry Smith 
1492758f92a0SBarry Smith    Output Parameters:
1493758f92a0SBarry Smith .  a   - array to hold history
1494758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1495758f92a0SBarry Smith          negative if not converged) for each solve.
1496758f92a0SBarry Smith -  na  - size of a and its
1497758f92a0SBarry Smith 
1498758f92a0SBarry Smith    Notes:
1499758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1500758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1501758f92a0SBarry Smith 
1502758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1503758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1504758f92a0SBarry Smith    during the section of code that is being timed.
1505758f92a0SBarry Smith 
1506758f92a0SBarry Smith    Level: intermediate
1507758f92a0SBarry Smith 
1508758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1509758f92a0SBarry Smith 
1510758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1511758f92a0SBarry Smith 
1512758f92a0SBarry Smith @*/
151363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1514758f92a0SBarry Smith {
1515758f92a0SBarry Smith   PetscFunctionBegin;
15164482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1517758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1518758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1519758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
15203a40ed3dSBarry Smith   PetscFunctionReturn(0);
1521c9005455SLois Curfman McInnes }
1522c9005455SLois Curfman McInnes 
1523e74ef692SMatthew Knepley #undef __FUNCT__
1524e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1525ac226902SBarry Smith /*@C
152676b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
152776b2cf59SMatthew Knepley   at the beginning of every step of the iteration.
152876b2cf59SMatthew Knepley 
152976b2cf59SMatthew Knepley   Collective on SNES
153076b2cf59SMatthew Knepley 
153176b2cf59SMatthew Knepley   Input Parameters:
153276b2cf59SMatthew Knepley . snes - The nonlinear solver context
153376b2cf59SMatthew Knepley . func - The function
153476b2cf59SMatthew Knepley 
153576b2cf59SMatthew Knepley   Calling sequence of func:
1536b5d30489SBarry Smith . func (SNES snes, PetscInt step);
153776b2cf59SMatthew Knepley 
153876b2cf59SMatthew Knepley . step - The current step of the iteration
153976b2cf59SMatthew Knepley 
154076b2cf59SMatthew Knepley   Level: intermediate
154176b2cf59SMatthew Knepley 
154276b2cf59SMatthew Knepley .keywords: SNES, update
1543b5d30489SBarry Smith 
154476b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
154576b2cf59SMatthew Knepley @*/
154663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
154776b2cf59SMatthew Knepley {
154876b2cf59SMatthew Knepley   PetscFunctionBegin;
15494482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
155076b2cf59SMatthew Knepley   snes->update = func;
155176b2cf59SMatthew Knepley   PetscFunctionReturn(0);
155276b2cf59SMatthew Knepley }
155376b2cf59SMatthew Knepley 
1554e74ef692SMatthew Knepley #undef __FUNCT__
1555e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
155676b2cf59SMatthew Knepley /*@
155776b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
155876b2cf59SMatthew Knepley 
155976b2cf59SMatthew Knepley   Not collective
156076b2cf59SMatthew Knepley 
156176b2cf59SMatthew Knepley   Input Parameters:
156276b2cf59SMatthew Knepley . snes - The nonlinear solver context
156376b2cf59SMatthew Knepley . step - The current step of the iteration
156476b2cf59SMatthew Knepley 
1565205452f4SMatthew Knepley   Level: intermediate
1566205452f4SMatthew Knepley 
156776b2cf59SMatthew Knepley .keywords: SNES, update
156876b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
156976b2cf59SMatthew Knepley @*/
157063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
157176b2cf59SMatthew Knepley {
157276b2cf59SMatthew Knepley   PetscFunctionBegin;
157376b2cf59SMatthew Knepley   PetscFunctionReturn(0);
157476b2cf59SMatthew Knepley }
157576b2cf59SMatthew Knepley 
15764a2ae208SSatish Balay #undef __FUNCT__
15774a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
15789b94acceSBarry Smith /*
15799b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15809b94acceSBarry Smith    positive parameter delta.
15819b94acceSBarry Smith 
15829b94acceSBarry Smith     Input Parameters:
1583c7afd0dbSLois Curfman McInnes +   snes - the SNES context
15849b94acceSBarry Smith .   y - approximate solution of linear system
15859b94acceSBarry Smith .   fnorm - 2-norm of current function
1586c7afd0dbSLois Curfman McInnes -   delta - trust region size
15879b94acceSBarry Smith 
15889b94acceSBarry Smith     Output Parameters:
1589c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
15909b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15919b94acceSBarry Smith     region, and exceeds zero otherwise.
1592c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
15939b94acceSBarry Smith 
15949b94acceSBarry Smith     Note:
15954b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
15969b94acceSBarry Smith     is set to be the maximum allowable step size.
15979b94acceSBarry Smith 
15989b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
15999b94acceSBarry Smith */
1600dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
16019b94acceSBarry Smith {
1602064f8208SBarry Smith   PetscReal      nrm;
1603ea709b57SSatish Balay   PetscScalar    cnorm;
1604dfbe8321SBarry Smith   PetscErrorCode ierr;
16053a40ed3dSBarry Smith 
16063a40ed3dSBarry Smith   PetscFunctionBegin;
16074482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16084482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1609c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1610184914b5SBarry Smith 
1611064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1612064f8208SBarry Smith   if (nrm > *delta) {
1613064f8208SBarry Smith      nrm = *delta/nrm;
1614064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1615064f8208SBarry Smith      cnorm = nrm;
16162dcb1b2aSMatthew Knepley      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
16179b94acceSBarry Smith      *ynorm = *delta;
16189b94acceSBarry Smith   } else {
16199b94acceSBarry Smith      *gpnorm = 0.0;
1620064f8208SBarry Smith      *ynorm = nrm;
16219b94acceSBarry Smith   }
16223a40ed3dSBarry Smith   PetscFunctionReturn(0);
16239b94acceSBarry Smith }
16249b94acceSBarry Smith 
16254a2ae208SSatish Balay #undef __FUNCT__
16264a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
16279b94acceSBarry Smith /*@
1628f69a0ea3SMatthew Knepley    SNESSolve - Solves a nonlinear system F(x) = b.
1629f69a0ea3SMatthew Knepley    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
16309b94acceSBarry Smith 
1631c7afd0dbSLois Curfman McInnes    Collective on SNES
1632c7afd0dbSLois Curfman McInnes 
1633b2002411SLois Curfman McInnes    Input Parameters:
1634c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1635f69a0ea3SMatthew Knepley .  b - the constant part of the equation, or PETSC_NULL to use zero.
163670e92668SMatthew Knepley -  x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution()
16379b94acceSBarry Smith 
1638b2002411SLois Curfman McInnes    Notes:
16398ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
16408ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16418ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16428ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16438ddd3da0SLois Curfman McInnes 
164436851e7fSLois Curfman McInnes    Level: beginner
164536851e7fSLois Curfman McInnes 
16469b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16479b94acceSBarry Smith 
164870e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution()
16499b94acceSBarry Smith @*/
1650f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
16519b94acceSBarry Smith {
1652dfbe8321SBarry Smith   PetscErrorCode ierr;
1653f1af5d2fSBarry Smith   PetscTruth     flg;
1654052efed2SBarry Smith 
16553a40ed3dSBarry Smith   PetscFunctionBegin;
16564482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16571302d50aSBarry Smith   if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
165874637425SBarry Smith 
1659f69a0ea3SMatthew Knepley   if (b) {
1660f69a0ea3SMatthew Knepley     ierr = SNESSetRhs(snes, b); CHKERRQ(ierr);
16611096aae1SMatthew Knepley     if (!snes->vec_func) {
16621096aae1SMatthew Knepley       Vec r;
16631096aae1SMatthew Knepley 
16641096aae1SMatthew Knepley       ierr = VecDuplicate(b, &r); CHKERRQ(ierr);
16651096aae1SMatthew Knepley       ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr);
16661096aae1SMatthew Knepley     }
1667f69a0ea3SMatthew Knepley   }
166870e92668SMatthew Knepley   if (x) {
1669f69a0ea3SMatthew Knepley     PetscValidHeaderSpecific(x,VEC_COOKIE,3);
1670f69a0ea3SMatthew Knepley     PetscCheckSameComm(snes,1,x,3);
167170e92668SMatthew Knepley   } else {
167270e92668SMatthew Knepley     ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr);
167370e92668SMatthew Knepley     if (!x) {
167470e92668SMatthew Knepley       ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr);
167570e92668SMatthew Knepley     }
167670e92668SMatthew Knepley   }
167770e92668SMatthew Knepley   snes->vec_sol = snes->vec_sol_always = x;
167870e92668SMatthew Knepley   if (!snes->setupcalled) {
167970e92668SMatthew Knepley     ierr = SNESSetUp(snes);CHKERRQ(ierr);
168070e92668SMatthew Knepley   }
1681abc0a331SBarry Smith   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
1682d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
168350ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1684d5e45103SBarry Smith 
1685d5e45103SBarry Smith   ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN);
1686d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
168738f152feSBarry Smith     /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */
168819717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr);
1689d5e45103SBarry Smith   } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) {
1690d5e45103SBarry Smith     /* translate exception into SNES not converged reason */
1691d5e45103SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
169238f152feSBarry Smith     ierr = 0;
169338f152feSBarry Smith   }
169438f152feSBarry Smith   CHKERRQ(ierr);
1695d5e45103SBarry Smith 
1696d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1697b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
16988b179ff8SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); }
1699da9b6338SBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1700da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
17015968eb51SBarry Smith   if (snes->printreason) {
17025968eb51SBarry Smith     if (snes->reason > 0) {
17039dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17045968eb51SBarry Smith     } else {
17059dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17065968eb51SBarry Smith     }
17075968eb51SBarry Smith   }
17085968eb51SBarry Smith 
17093a40ed3dSBarry Smith   PetscFunctionReturn(0);
17109b94acceSBarry Smith }
17119b94acceSBarry Smith 
17129b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
17139b94acceSBarry Smith 
17144a2ae208SSatish Balay #undef __FUNCT__
17154a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
171682bf6240SBarry Smith /*@C
17174b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
17189b94acceSBarry Smith 
1719fee21e36SBarry Smith    Collective on SNES
1720fee21e36SBarry Smith 
1721c7afd0dbSLois Curfman McInnes    Input Parameters:
1722c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1723454a90a3SBarry Smith -  type - a known method
1724c7afd0dbSLois Curfman McInnes 
1725c7afd0dbSLois Curfman McInnes    Options Database Key:
1726454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1727c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1728ae12b187SLois Curfman McInnes 
17299b94acceSBarry Smith    Notes:
1730e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
17314b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
1732c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17334b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
1734c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17359b94acceSBarry Smith 
1736ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1737ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1738ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1739ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1740ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1741ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1742ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1743ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1744ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1745b0a32e0cSBarry Smith   appropriate method.
174636851e7fSLois Curfman McInnes 
174736851e7fSLois Curfman McInnes   Level: intermediate
1748a703fe33SLois Curfman McInnes 
1749454a90a3SBarry Smith .keywords: SNES, set, type
1750435da068SBarry Smith 
1751435da068SBarry Smith .seealso: SNESType, SNESCreate()
1752435da068SBarry Smith 
17539b94acceSBarry Smith @*/
1754e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type)
17559b94acceSBarry Smith {
1756dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
17576831982aSBarry Smith   PetscTruth     match;
17583a40ed3dSBarry Smith 
17593a40ed3dSBarry Smith   PetscFunctionBegin;
17604482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17614482741eSBarry Smith   PetscValidCharPointer(type,2);
176282bf6240SBarry Smith 
17636831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
17640f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
176582bf6240SBarry Smith 
176682bf6240SBarry Smith   if (snes->setupcalled) {
17674dc4c822SBarry Smith     snes->setupcalled = PETSC_FALSE;
1768e1311b90SBarry Smith     ierr              = (*(snes)->destroy)(snes);CHKERRQ(ierr);
176982bf6240SBarry Smith     snes->data        = 0;
177082bf6240SBarry Smith   }
17717d1a2b2bSBarry Smith 
17729b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
177382bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1774b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1775958c9bccSBarry Smith   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
1776606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
177782bf6240SBarry Smith   snes->data = 0;
17783a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1779454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
17803a40ed3dSBarry Smith   PetscFunctionReturn(0);
17819b94acceSBarry Smith }
17829b94acceSBarry Smith 
1783a847f771SSatish Balay 
17849b94acceSBarry Smith /* --------------------------------------------------------------------- */
17854a2ae208SSatish Balay #undef __FUNCT__
17864a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
17879b94acceSBarry Smith /*@C
17889b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1789f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
17909b94acceSBarry Smith 
1791fee21e36SBarry Smith    Not Collective
1792fee21e36SBarry Smith 
179336851e7fSLois Curfman McInnes    Level: advanced
179436851e7fSLois Curfman McInnes 
17959b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
17969b94acceSBarry Smith 
17979b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
17989b94acceSBarry Smith @*/
179963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
18009b94acceSBarry Smith {
1801dfbe8321SBarry Smith   PetscErrorCode ierr;
180282bf6240SBarry Smith 
18033a40ed3dSBarry Smith   PetscFunctionBegin;
180482bf6240SBarry Smith   if (SNESList) {
1805b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
180682bf6240SBarry Smith     SNESList = 0;
18079b94acceSBarry Smith   }
18084c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
18093a40ed3dSBarry Smith   PetscFunctionReturn(0);
18109b94acceSBarry Smith }
18119b94acceSBarry Smith 
18124a2ae208SSatish Balay #undef __FUNCT__
18134a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
18149b94acceSBarry Smith /*@C
18159a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
18169b94acceSBarry Smith 
1817c7afd0dbSLois Curfman McInnes    Not Collective
1818c7afd0dbSLois Curfman McInnes 
18199b94acceSBarry Smith    Input Parameter:
18204b0e389bSBarry Smith .  snes - nonlinear solver context
18219b94acceSBarry Smith 
18229b94acceSBarry Smith    Output Parameter:
18233a7fca6bSBarry Smith .  type - SNES method (a character string)
18249b94acceSBarry Smith 
182536851e7fSLois Curfman McInnes    Level: intermediate
182636851e7fSLois Curfman McInnes 
1827454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
18289b94acceSBarry Smith @*/
182963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type)
18309b94acceSBarry Smith {
18313a40ed3dSBarry Smith   PetscFunctionBegin;
18324482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18334482741eSBarry Smith   PetscValidPointer(type,2);
1834454a90a3SBarry Smith   *type = snes->type_name;
18353a40ed3dSBarry Smith   PetscFunctionReturn(0);
18369b94acceSBarry Smith }
18379b94acceSBarry Smith 
18384a2ae208SSatish Balay #undef __FUNCT__
18394a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
18409b94acceSBarry Smith /*@C
18419b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
18429b94acceSBarry Smith    stored.
18439b94acceSBarry Smith 
1844c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1845c7afd0dbSLois Curfman McInnes 
18469b94acceSBarry Smith    Input Parameter:
18479b94acceSBarry Smith .  snes - the SNES context
18489b94acceSBarry Smith 
18499b94acceSBarry Smith    Output Parameter:
18509b94acceSBarry Smith .  x - the solution
18519b94acceSBarry Smith 
185270e92668SMatthew Knepley    Level: intermediate
185336851e7fSLois Curfman McInnes 
18549b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18559b94acceSBarry Smith 
185670e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
18579b94acceSBarry Smith @*/
185863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
18599b94acceSBarry Smith {
18603a40ed3dSBarry Smith   PetscFunctionBegin;
18614482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18624482741eSBarry Smith   PetscValidPointer(x,2);
18639b94acceSBarry Smith   *x = snes->vec_sol_always;
18643a40ed3dSBarry Smith   PetscFunctionReturn(0);
18659b94acceSBarry Smith }
18669b94acceSBarry Smith 
18674a2ae208SSatish Balay #undef __FUNCT__
186870e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution"
186970e92668SMatthew Knepley /*@C
187070e92668SMatthew Knepley    SNESSetSolution - Sets the vector where the approximate solution is stored.
187170e92668SMatthew Knepley 
187270e92668SMatthew Knepley    Not Collective, but Vec is parallel if SNES is parallel
187370e92668SMatthew Knepley 
187470e92668SMatthew Knepley    Input Parameters:
187570e92668SMatthew Knepley +  snes - the SNES context
187670e92668SMatthew Knepley -  x - the solution
187770e92668SMatthew Knepley 
187870e92668SMatthew Knepley    Output Parameter:
187970e92668SMatthew Knepley 
188070e92668SMatthew Knepley    Level: intermediate
188170e92668SMatthew Knepley 
188270e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution
188370e92668SMatthew Knepley 
188470e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
188570e92668SMatthew Knepley @*/
188670e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x)
188770e92668SMatthew Knepley {
188870e92668SMatthew Knepley   PetscFunctionBegin;
188970e92668SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
189070e92668SMatthew Knepley   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
189170e92668SMatthew Knepley   PetscCheckSameComm(snes,1,x,2);
189270e92668SMatthew Knepley   snes->vec_sol_always = x;
189370e92668SMatthew Knepley   PetscFunctionReturn(0);
189470e92668SMatthew Knepley }
189570e92668SMatthew Knepley 
189670e92668SMatthew Knepley #undef __FUNCT__
18974a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
18989b94acceSBarry Smith /*@C
18999b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
19009b94acceSBarry Smith    stored.
19019b94acceSBarry Smith 
1902c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1903c7afd0dbSLois Curfman McInnes 
19049b94acceSBarry Smith    Input Parameter:
19059b94acceSBarry Smith .  snes - the SNES context
19069b94acceSBarry Smith 
19079b94acceSBarry Smith    Output Parameter:
19089b94acceSBarry Smith .  x - the solution update
19099b94acceSBarry Smith 
191036851e7fSLois Curfman McInnes    Level: advanced
191136851e7fSLois Curfman McInnes 
19129b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
19139b94acceSBarry Smith 
19149b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
19159b94acceSBarry Smith @*/
191663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
19179b94acceSBarry Smith {
19183a40ed3dSBarry Smith   PetscFunctionBegin;
19194482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19204482741eSBarry Smith   PetscValidPointer(x,2);
19219b94acceSBarry Smith   *x = snes->vec_sol_update_always;
19223a40ed3dSBarry Smith   PetscFunctionReturn(0);
19239b94acceSBarry Smith }
19249b94acceSBarry Smith 
19254a2ae208SSatish Balay #undef __FUNCT__
19264a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
19279b94acceSBarry Smith /*@C
19283638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19299b94acceSBarry Smith 
1930c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1931c7afd0dbSLois Curfman McInnes 
19329b94acceSBarry Smith    Input Parameter:
19339b94acceSBarry Smith .  snes - the SNES context
19349b94acceSBarry Smith 
19359b94acceSBarry Smith    Output Parameter:
19367bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
193770e92668SMatthew Knepley .  func - the function (or PETSC_NULL)
193870e92668SMatthew Knepley -  ctx - the function context (or PETSC_NULL)
19399b94acceSBarry Smith 
194036851e7fSLois Curfman McInnes    Level: advanced
194136851e7fSLois Curfman McInnes 
1942a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19439b94acceSBarry Smith 
19444b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
19459b94acceSBarry Smith @*/
194670e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
19479b94acceSBarry Smith {
19483a40ed3dSBarry Smith   PetscFunctionBegin;
19494482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19507bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
195100036973SBarry Smith   if (func) *func = snes->computefunction;
195270e92668SMatthew Knepley   if (ctx)  *ctx  = snes->funP;
19533a40ed3dSBarry Smith   PetscFunctionReturn(0);
19549b94acceSBarry Smith }
19559b94acceSBarry Smith 
19564a2ae208SSatish Balay #undef __FUNCT__
19574a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
19583c7409f5SSatish Balay /*@C
19593c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1960d850072dSLois Curfman McInnes    SNES options in the database.
19613c7409f5SSatish Balay 
1962fee21e36SBarry Smith    Collective on SNES
1963fee21e36SBarry Smith 
1964c7afd0dbSLois Curfman McInnes    Input Parameter:
1965c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1966c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
1967c7afd0dbSLois Curfman McInnes 
1968d850072dSLois Curfman McInnes    Notes:
1969a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1970c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
1971d850072dSLois Curfman McInnes 
197236851e7fSLois Curfman McInnes    Level: advanced
197336851e7fSLois Curfman McInnes 
19743c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1975a86d99e1SLois Curfman McInnes 
1976a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
19773c7409f5SSatish Balay @*/
197863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
19793c7409f5SSatish Balay {
1980dfbe8321SBarry Smith   PetscErrorCode ierr;
19813c7409f5SSatish Balay 
19823a40ed3dSBarry Smith   PetscFunctionBegin;
19834482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1984639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
198594b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
19863a40ed3dSBarry Smith   PetscFunctionReturn(0);
19873c7409f5SSatish Balay }
19883c7409f5SSatish Balay 
19894a2ae208SSatish Balay #undef __FUNCT__
19904a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
19913c7409f5SSatish Balay /*@C
1992f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1993d850072dSLois Curfman McInnes    SNES options in the database.
19943c7409f5SSatish Balay 
1995fee21e36SBarry Smith    Collective on SNES
1996fee21e36SBarry Smith 
1997c7afd0dbSLois Curfman McInnes    Input Parameters:
1998c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1999c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2000c7afd0dbSLois Curfman McInnes 
2001d850072dSLois Curfman McInnes    Notes:
2002a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2003c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2004d850072dSLois Curfman McInnes 
200536851e7fSLois Curfman McInnes    Level: advanced
200636851e7fSLois Curfman McInnes 
20073c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2008a86d99e1SLois Curfman McInnes 
2009a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20103c7409f5SSatish Balay @*/
201163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
20123c7409f5SSatish Balay {
2013dfbe8321SBarry Smith   PetscErrorCode ierr;
20143c7409f5SSatish Balay 
20153a40ed3dSBarry Smith   PetscFunctionBegin;
20164482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2017639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
201894b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20193a40ed3dSBarry Smith   PetscFunctionReturn(0);
20203c7409f5SSatish Balay }
20213c7409f5SSatish Balay 
20224a2ae208SSatish Balay #undef __FUNCT__
20234a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
20249ab63eb5SSatish Balay /*@C
20253c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20263c7409f5SSatish Balay    SNES options in the database.
20273c7409f5SSatish Balay 
2028c7afd0dbSLois Curfman McInnes    Not Collective
2029c7afd0dbSLois Curfman McInnes 
20303c7409f5SSatish Balay    Input Parameter:
20313c7409f5SSatish Balay .  snes - the SNES context
20323c7409f5SSatish Balay 
20333c7409f5SSatish Balay    Output Parameter:
20343c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20353c7409f5SSatish Balay 
20369ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
20379ab63eb5SSatish Balay    sufficient length to hold the prefix.
20389ab63eb5SSatish Balay 
203936851e7fSLois Curfman McInnes    Level: advanced
204036851e7fSLois Curfman McInnes 
20413c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2042a86d99e1SLois Curfman McInnes 
2043a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20443c7409f5SSatish Balay @*/
2045e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
20463c7409f5SSatish Balay {
2047dfbe8321SBarry Smith   PetscErrorCode ierr;
20483c7409f5SSatish Balay 
20493a40ed3dSBarry Smith   PetscFunctionBegin;
20504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2051639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
20523a40ed3dSBarry Smith   PetscFunctionReturn(0);
20533c7409f5SSatish Balay }
20543c7409f5SSatish Balay 
2055b2002411SLois Curfman McInnes 
20564a2ae208SSatish Balay #undef __FUNCT__
20574a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
20583cea93caSBarry Smith /*@C
20593cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
20603cea93caSBarry Smith 
20617f6c08e0SMatthew Knepley   Level: advanced
20623cea93caSBarry Smith @*/
206363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2064b2002411SLois Curfman McInnes {
2065e2d1d2b7SBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
2066dfbe8321SBarry Smith   PetscErrorCode ierr;
2067b2002411SLois Curfman McInnes 
2068b2002411SLois Curfman McInnes   PetscFunctionBegin;
2069b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2070c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2071b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2072b2002411SLois Curfman McInnes }
2073da9b6338SBarry Smith 
2074da9b6338SBarry Smith #undef __FUNCT__
2075da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
207663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2077da9b6338SBarry Smith {
2078dfbe8321SBarry Smith   PetscErrorCode ierr;
207977431f27SBarry Smith   PetscInt       N,i,j;
2080da9b6338SBarry Smith   Vec            u,uh,fh;
2081da9b6338SBarry Smith   PetscScalar    value;
2082da9b6338SBarry Smith   PetscReal      norm;
2083da9b6338SBarry Smith 
2084da9b6338SBarry Smith   PetscFunctionBegin;
2085da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2086da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2087da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2088da9b6338SBarry Smith 
2089da9b6338SBarry Smith   /* currently only works for sequential */
2090da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2091da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2092da9b6338SBarry Smith   for (i=0; i<N; i++) {
2093da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
209477431f27SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2095da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2096ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2097da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
20983ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2099da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
210077431f27SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2101da9b6338SBarry Smith       value = -value;
2102da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2103da9b6338SBarry Smith     }
2104da9b6338SBarry Smith   }
2105da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2106da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2107da9b6338SBarry Smith   PetscFunctionReturn(0);
2108da9b6338SBarry Smith }
2109