xref: /petsc/src/snes/interface/snes.c (revision dc67937b7f06f2e4bb3627d61614c9b821dfc2de)
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"
51452baeb72SSatish Balay /*@
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"
55652baeb72SSatish Balay /*@
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");
798e9a2bbcdSBarry Smith     CHKMEMQ;
79919717074SBarry Smith     ierr = (*snes->computefunction)(snes,x,y,snes->funP);
800e9a2bbcdSBarry 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) {
812016dedfbSBarry Smith     ierr = VecAXPY(y,-1.0,snes->afine);CHKERRQ(ierr);
8133ab0aad5SBarry Smith   }
814ae3c334cSLois Curfman McInnes   snes->nfuncs++;
815d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8163a40ed3dSBarry Smith   PetscFunctionReturn(0);
8179b94acceSBarry Smith }
8189b94acceSBarry Smith 
8194a2ae208SSatish Balay #undef __FUNCT__
8204a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
82162fef451SLois Curfman McInnes /*@
82262fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
82362fef451SLois Curfman McInnes    set with SNESSetJacobian().
82462fef451SLois Curfman McInnes 
825c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
826c7afd0dbSLois Curfman McInnes 
82762fef451SLois Curfman McInnes    Input Parameters:
828c7afd0dbSLois Curfman McInnes +  snes - the SNES context
829c7afd0dbSLois Curfman McInnes -  x - input vector
83062fef451SLois Curfman McInnes 
83162fef451SLois Curfman McInnes    Output Parameters:
832c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
83362fef451SLois Curfman McInnes .  B - optional preconditioning matrix
834c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
835fee21e36SBarry Smith 
83662fef451SLois Curfman McInnes    Notes:
83762fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
83862fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
83962fef451SLois Curfman McInnes 
84094b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
841dc5a77f8SLois Curfman McInnes    flag parameter.
84262fef451SLois Curfman McInnes 
84336851e7fSLois Curfman McInnes    Level: developer
84436851e7fSLois Curfman McInnes 
84562fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
84662fef451SLois Curfman McInnes 
84794b7f48cSBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators()
84862fef451SLois Curfman McInnes @*/
84963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8509b94acceSBarry Smith {
851dfbe8321SBarry Smith   PetscErrorCode ierr;
8523a40ed3dSBarry Smith 
8533a40ed3dSBarry Smith   PetscFunctionBegin;
8544482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8554482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
8564482741eSBarry Smith   PetscValidPointer(flg,5);
857c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
8583a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
859d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
860c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
861d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
862*dc67937bSBarry Smith   CHKMEMQ;
8639b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
864*dc67937bSBarry Smith   CHKMEMQ;
865d64ed03dSBarry Smith   PetscStackPop;
866d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
8676d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
8684482741eSBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
8694482741eSBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE,4);
8703a40ed3dSBarry Smith   PetscFunctionReturn(0);
8719b94acceSBarry Smith }
8729b94acceSBarry Smith 
8734a2ae208SSatish Balay #undef __FUNCT__
8744a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
8759b94acceSBarry Smith /*@C
8769b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
877044dda88SLois Curfman McInnes    location to store the matrix.
8789b94acceSBarry Smith 
879c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
880c7afd0dbSLois Curfman McInnes 
8819b94acceSBarry Smith    Input Parameters:
882c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8839b94acceSBarry Smith .  A - Jacobian matrix
8849b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8859b94acceSBarry Smith .  func - Jacobian evaluation routine
886c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
8872cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
8889b94acceSBarry Smith 
8899b94acceSBarry Smith    Calling sequence of func:
8908d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
8919b94acceSBarry Smith 
892c7afd0dbSLois Curfman McInnes +  x - input vector
8939b94acceSBarry Smith .  A - Jacobian matrix
8949b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
895ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
89694b7f48cSBarry Smith    structure (same as flag in KSPSetOperators())
897c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
8989b94acceSBarry Smith 
8999b94acceSBarry Smith    Notes:
90094b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
9012cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
902ac21db08SLois Curfman McInnes 
903ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
9049b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
9059b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
9069b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
9079b94acceSBarry Smith    throughout the global iterations.
9089b94acceSBarry Smith 
90936851e7fSLois Curfman McInnes    Level: beginner
91036851e7fSLois Curfman McInnes 
9119b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
9129b94acceSBarry Smith 
9133960baedSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor()
9149b94acceSBarry Smith @*/
91563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
9169b94acceSBarry Smith {
917dfbe8321SBarry Smith   PetscErrorCode ierr;
9183a7fca6bSBarry Smith 
9193a40ed3dSBarry Smith   PetscFunctionBegin;
9204482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9214482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
9224482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
923c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
924c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
9253a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
9263a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
9273a7fca6bSBarry Smith   if (A) {
9283a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
9299b94acceSBarry Smith     snes->jacobian = A;
9303a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
9313a7fca6bSBarry Smith   }
9323a7fca6bSBarry Smith   if (B) {
9333a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
9349b94acceSBarry Smith     snes->jacobian_pre = B;
9353a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
9363a7fca6bSBarry Smith   }
93769a612a9SBarry Smith   ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
9383a40ed3dSBarry Smith   PetscFunctionReturn(0);
9399b94acceSBarry Smith }
94062fef451SLois Curfman McInnes 
9414a2ae208SSatish Balay #undef __FUNCT__
9424a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
943c2aafc4cSSatish Balay /*@C
944b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
945b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
946b4fd4287SBarry Smith 
947c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
948c7afd0dbSLois Curfman McInnes 
949b4fd4287SBarry Smith    Input Parameter:
950b4fd4287SBarry Smith .  snes - the nonlinear solver context
951b4fd4287SBarry Smith 
952b4fd4287SBarry Smith    Output Parameters:
953c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
954b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
95570e92668SMatthew Knepley .  func - location to put Jacobian function (or PETSC_NULL)
95670e92668SMatthew Knepley -  ctx - location to stash Jacobian ctx (or PETSC_NULL)
957fee21e36SBarry Smith 
95836851e7fSLois Curfman McInnes    Level: advanced
95936851e7fSLois Curfman McInnes 
960b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
961b4fd4287SBarry Smith @*/
96270e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx)
963b4fd4287SBarry Smith {
9643a40ed3dSBarry Smith   PetscFunctionBegin;
9654482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
966b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
967b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
96800036973SBarry Smith   if (func) *func = snes->computejacobian;
96970e92668SMatthew Knepley   if (ctx)  *ctx  = snes->jacP;
9703a40ed3dSBarry Smith   PetscFunctionReturn(0);
971b4fd4287SBarry Smith }
972b4fd4287SBarry Smith 
9739b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
97463dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
9759b94acceSBarry Smith 
9764a2ae208SSatish Balay #undef __FUNCT__
9774a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
9789b94acceSBarry Smith /*@
9799b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
980272ac6f2SLois Curfman McInnes    of a nonlinear solver.
9819b94acceSBarry Smith 
982fee21e36SBarry Smith    Collective on SNES
983fee21e36SBarry Smith 
984c7afd0dbSLois Curfman McInnes    Input Parameters:
98570e92668SMatthew Knepley .  snes - the SNES context
986c7afd0dbSLois Curfman McInnes 
987272ac6f2SLois Curfman McInnes    Notes:
988272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
989272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
990272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
991272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
992272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
993272ac6f2SLois Curfman McInnes 
99436851e7fSLois Curfman McInnes    Level: advanced
99536851e7fSLois Curfman McInnes 
9969b94acceSBarry Smith .keywords: SNES, nonlinear, setup
9979b94acceSBarry Smith 
9989b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
9999b94acceSBarry Smith @*/
100070e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes)
10019b94acceSBarry Smith {
1002dfbe8321SBarry Smith   PetscErrorCode ierr;
10034b27c08aSLois Curfman McInnes   PetscTruth     flg, iseqtr;
10043a40ed3dSBarry Smith 
10053a40ed3dSBarry Smith   PetscFunctionBegin;
10064482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10074dc4c822SBarry Smith   if (snes->setupcalled) PetscFunctionReturn(0);
10089b94acceSBarry Smith 
1009b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1010c1f60f51SBarry Smith   /*
1011c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1012dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1013c1f60f51SBarry Smith   */
1014c1f60f51SBarry Smith   if (flg) {
1015c1f60f51SBarry Smith     Mat J;
10165a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10175a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
101863ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr);
10193a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
10203a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1021c1f60f51SBarry Smith   }
102245fc7adcSBarry Smith 
1023c8d321e0SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE)
102445fc7adcSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
102545fc7adcSBarry Smith   if (flg) {
102645fc7adcSBarry Smith     Mat J;
102745fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
102845fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
102945fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
103045fc7adcSBarry Smith   }
103132a4b47aSMatthew Knepley #endif
103245fc7adcSBarry Smith 
1033b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1034c1f60f51SBarry Smith   /*
1035dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1036c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1037c1f60f51SBarry Smith    */
103831615d04SLois Curfman McInnes   if (flg) {
1039272ac6f2SLois Curfman McInnes     Mat  J;
1040b5d62d44SBarry Smith     KSP ksp;
104194b7f48cSBarry Smith     PC   pc;
1042f3ef73ceSBarry Smith 
10435a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10443a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
104563ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr);
10463a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
10473a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
10483a7fca6bSBarry Smith 
1049f3ef73ceSBarry Smith     /* force no preconditioner */
105094b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1051b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1052a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1053a9815358SBarry Smith     if (!flg) {
1054f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1055272ac6f2SLois Curfman McInnes     }
1056a9815358SBarry Smith   }
1057f3ef73ceSBarry Smith 
10581096aae1SMatthew Knepley   if (!snes->vec_func && !snes->afine) {
10591096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10601096aae1SMatthew Knepley   }
10611096aae1SMatthew Knepley   if (!snes->computefunction && !snes->afine) {
10621096aae1SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
10631096aae1SMatthew Knepley   }
106429bbc08cSBarry Smith   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1065a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
106629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1067a8c6a408SBarry Smith   }
1068a703fe33SLois Curfman McInnes 
1069a703fe33SLois Curfman McInnes   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
10704b27c08aSLois Curfman McInnes   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
10716831982aSBarry Smith   if (snes->ksp_ewconv && !iseqtr) {
107294b7f48cSBarry Smith     KSP ksp;
107394b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1074186905e3SBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1075a703fe33SLois Curfman McInnes   }
10764b27c08aSLois Curfman McInnes 
1077a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
10787aaed0d8SBarry Smith   snes->setupcalled = PETSC_TRUE;
10793a40ed3dSBarry Smith   PetscFunctionReturn(0);
10809b94acceSBarry Smith }
10819b94acceSBarry Smith 
10824a2ae208SSatish Balay #undef __FUNCT__
10834a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
108452baeb72SSatish Balay /*@
10859b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
10869b94acceSBarry Smith    with SNESCreate().
10879b94acceSBarry Smith 
1088c7afd0dbSLois Curfman McInnes    Collective on SNES
1089c7afd0dbSLois Curfman McInnes 
10909b94acceSBarry Smith    Input Parameter:
10919b94acceSBarry Smith .  snes - the SNES context
10929b94acceSBarry Smith 
109336851e7fSLois Curfman McInnes    Level: beginner
109436851e7fSLois Curfman McInnes 
10959b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
10969b94acceSBarry Smith 
109763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
10989b94acceSBarry Smith @*/
109963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes)
11009b94acceSBarry Smith {
110177431f27SBarry Smith   PetscInt       i;
11026849ba73SBarry Smith   PetscErrorCode ierr;
11033a40ed3dSBarry Smith 
11043a40ed3dSBarry Smith   PetscFunctionBegin;
11054482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
11063a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1107d4bb536fSBarry Smith 
1108be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
11090f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1110be0abb6dSBarry Smith 
1111e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1112606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
11133a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
11143a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
11153ab0aad5SBarry Smith   if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);}
111694b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1117522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1118b8a78c4aSBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1119b8a78c4aSBarry Smith     if (snes->monitordestroy[i]) {
1120b8a78c4aSBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1121b8a78c4aSBarry Smith     }
1122b8a78c4aSBarry Smith   }
1123a79aaaedSSatish Balay   ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr);
11243a40ed3dSBarry Smith   PetscFunctionReturn(0);
11259b94acceSBarry Smith }
11269b94acceSBarry Smith 
11279b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11289b94acceSBarry Smith 
11294a2ae208SSatish Balay #undef __FUNCT__
11304a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
11319b94acceSBarry Smith /*@
1132d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11339b94acceSBarry Smith 
1134c7afd0dbSLois Curfman McInnes    Collective on SNES
1135c7afd0dbSLois Curfman McInnes 
11369b94acceSBarry Smith    Input Parameters:
1137c7afd0dbSLois Curfman McInnes +  snes - the SNES context
113870441072SBarry Smith .  abstol - absolute convergence tolerance
113933174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
114033174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
114133174efeSLois Curfman McInnes            of the change in the solution between steps
114233174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1143c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1144fee21e36SBarry Smith 
114533174efeSLois Curfman McInnes    Options Database Keys:
114670441072SBarry Smith +    -snes_atol <abstol> - Sets abstol
1147c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1148c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1149c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1150c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
11519b94acceSBarry Smith 
1152d7a720efSLois Curfman McInnes    Notes:
11539b94acceSBarry Smith    The default maximum number of iterations is 50.
11549b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11559b94acceSBarry Smith 
115636851e7fSLois Curfman McInnes    Level: intermediate
115736851e7fSLois Curfman McInnes 
115833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11599b94acceSBarry Smith 
11602492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance()
11619b94acceSBarry Smith @*/
116263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf)
11639b94acceSBarry Smith {
11643a40ed3dSBarry Smith   PetscFunctionBegin;
11654482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
116670441072SBarry Smith   if (abstol != PETSC_DEFAULT)  snes->abstol      = abstol;
1167d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1168d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1169d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1170d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11713a40ed3dSBarry Smith   PetscFunctionReturn(0);
11729b94acceSBarry Smith }
11739b94acceSBarry Smith 
11744a2ae208SSatish Balay #undef __FUNCT__
11754a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
11769b94acceSBarry Smith /*@
117733174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
117833174efeSLois Curfman McInnes 
1179c7afd0dbSLois Curfman McInnes    Not Collective
1180c7afd0dbSLois Curfman McInnes 
118133174efeSLois Curfman McInnes    Input Parameters:
1182c7afd0dbSLois Curfman McInnes +  snes - the SNES context
118370441072SBarry Smith .  abstol - absolute convergence tolerance
118433174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
118533174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
118633174efeSLois Curfman McInnes            of the change in the solution between steps
118733174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1188c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1189fee21e36SBarry Smith 
119033174efeSLois Curfman McInnes    Notes:
119133174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
119233174efeSLois Curfman McInnes 
119336851e7fSLois Curfman McInnes    Level: intermediate
119436851e7fSLois Curfman McInnes 
119533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
119633174efeSLois Curfman McInnes 
119733174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
119833174efeSLois Curfman McInnes @*/
119963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf)
120033174efeSLois Curfman McInnes {
12013a40ed3dSBarry Smith   PetscFunctionBegin;
12024482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
120370441072SBarry Smith   if (abstol)  *abstol  = snes->abstol;
120433174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
120533174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
120633174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
120733174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
12083a40ed3dSBarry Smith   PetscFunctionReturn(0);
120933174efeSLois Curfman McInnes }
121033174efeSLois Curfman McInnes 
12114a2ae208SSatish Balay #undef __FUNCT__
12124a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
121333174efeSLois Curfman McInnes /*@
12149b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
12159b94acceSBarry Smith 
1216fee21e36SBarry Smith    Collective on SNES
1217fee21e36SBarry Smith 
1218c7afd0dbSLois Curfman McInnes    Input Parameters:
1219c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1220c7afd0dbSLois Curfman McInnes -  tol - tolerance
1221c7afd0dbSLois Curfman McInnes 
12229b94acceSBarry Smith    Options Database Key:
1223c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
12249b94acceSBarry Smith 
122536851e7fSLois Curfman McInnes    Level: intermediate
122636851e7fSLois Curfman McInnes 
12279b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12289b94acceSBarry Smith 
12292492ecdbSBarry Smith .seealso: SNESSetTolerances()
12309b94acceSBarry Smith @*/
123163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
12329b94acceSBarry Smith {
12333a40ed3dSBarry Smith   PetscFunctionBegin;
12344482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12359b94acceSBarry Smith   snes->deltatol = tol;
12363a40ed3dSBarry Smith   PetscFunctionReturn(0);
12379b94acceSBarry Smith }
12389b94acceSBarry Smith 
1239df9fa365SBarry Smith /*
1240df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1241df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1242df9fa365SBarry Smith    macros instead of functions
1243df9fa365SBarry Smith */
12444a2ae208SSatish Balay #undef __FUNCT__
12454a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
124663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx)
1247ce1608b8SBarry Smith {
1248dfbe8321SBarry Smith   PetscErrorCode ierr;
1249ce1608b8SBarry Smith 
1250ce1608b8SBarry Smith   PetscFunctionBegin;
12514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1252ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1253ce1608b8SBarry Smith   PetscFunctionReturn(0);
1254ce1608b8SBarry Smith }
1255ce1608b8SBarry Smith 
12564a2ae208SSatish Balay #undef __FUNCT__
12574a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
125863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1259df9fa365SBarry Smith {
1260dfbe8321SBarry Smith   PetscErrorCode ierr;
1261df9fa365SBarry Smith 
1262df9fa365SBarry Smith   PetscFunctionBegin;
1263df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1264df9fa365SBarry Smith   PetscFunctionReturn(0);
1265df9fa365SBarry Smith }
1266df9fa365SBarry Smith 
12674a2ae208SSatish Balay #undef __FUNCT__
12684a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
126963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw)
1270df9fa365SBarry Smith {
1271dfbe8321SBarry Smith   PetscErrorCode ierr;
1272df9fa365SBarry Smith 
1273df9fa365SBarry Smith   PetscFunctionBegin;
1274df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1275df9fa365SBarry Smith   PetscFunctionReturn(0);
1276df9fa365SBarry Smith }
1277df9fa365SBarry Smith 
12789b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12799b94acceSBarry Smith 
12804a2ae208SSatish Balay #undef __FUNCT__
12814a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
12829b94acceSBarry Smith /*@C
12835cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
12849b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12859b94acceSBarry Smith    progress.
12869b94acceSBarry Smith 
1287fee21e36SBarry Smith    Collective on SNES
1288fee21e36SBarry Smith 
1289c7afd0dbSLois Curfman McInnes    Input Parameters:
1290c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1291c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1292b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1293b3006f0bSLois Curfman McInnes           monitor routine (use PETSC_NULL if no context is desitre)
1294b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1295b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
12969b94acceSBarry Smith 
1297c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1298a7cc72afSBarry Smith $     int func(SNES snes,PetscInt its, PetscReal norm,void *mctx)
1299c7afd0dbSLois Curfman McInnes 
1300c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1301c7afd0dbSLois Curfman McInnes .    its - iteration number
1302c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
130340a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
13049b94acceSBarry Smith 
13059665c990SLois Curfman McInnes    Options Database Keys:
1306c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1307c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1308c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1309c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1310c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1311c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1312c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1313c7afd0dbSLois Curfman McInnes                             the options database.
13149665c990SLois Curfman McInnes 
1315639f9d9dSBarry Smith    Notes:
13166bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
13176bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
13186bc08f3fSLois Curfman McInnes    order in which they were set.
1319639f9d9dSBarry Smith 
132036851e7fSLois Curfman McInnes    Level: intermediate
132136851e7fSLois Curfman McInnes 
13229b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
13239b94acceSBarry Smith 
13245cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
13259b94acceSBarry Smith @*/
132663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
13279b94acceSBarry Smith {
13283a40ed3dSBarry Smith   PetscFunctionBegin;
13294482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1330639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
133129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1332639f9d9dSBarry Smith   }
1333639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1334b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1335639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13363a40ed3dSBarry Smith   PetscFunctionReturn(0);
13379b94acceSBarry Smith }
13389b94acceSBarry Smith 
13394a2ae208SSatish Balay #undef __FUNCT__
13404a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
13415cd90555SBarry Smith /*@C
13425cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
13435cd90555SBarry Smith 
1344c7afd0dbSLois Curfman McInnes    Collective on SNES
1345c7afd0dbSLois Curfman McInnes 
13465cd90555SBarry Smith    Input Parameters:
13475cd90555SBarry Smith .  snes - the SNES context
13485cd90555SBarry Smith 
13491a480d89SAdministrator    Options Database Key:
1350c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1351c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1352c7afd0dbSLois Curfman McInnes     set via the options database
13535cd90555SBarry Smith 
13545cd90555SBarry Smith    Notes:
13555cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
13565cd90555SBarry Smith 
135736851e7fSLois Curfman McInnes    Level: intermediate
135836851e7fSLois Curfman McInnes 
13595cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
13605cd90555SBarry Smith 
13615cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
13625cd90555SBarry Smith @*/
136363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes)
13645cd90555SBarry Smith {
13655cd90555SBarry Smith   PetscFunctionBegin;
13664482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13675cd90555SBarry Smith   snes->numbermonitors = 0;
13685cd90555SBarry Smith   PetscFunctionReturn(0);
13695cd90555SBarry Smith }
13705cd90555SBarry Smith 
13714a2ae208SSatish Balay #undef __FUNCT__
13724a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
13739b94acceSBarry Smith /*@C
13749b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13759b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13769b94acceSBarry Smith 
1377fee21e36SBarry Smith    Collective on SNES
1378fee21e36SBarry Smith 
1379c7afd0dbSLois Curfman McInnes    Input Parameters:
1380c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1381c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1382c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1383c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
13849b94acceSBarry Smith 
1385c7afd0dbSLois Curfman McInnes    Calling sequence of func:
13866849ba73SBarry Smith $     PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1387c7afd0dbSLois Curfman McInnes 
1388c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1389c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1390184914b5SBarry Smith .    reason - reason for convergence/divergence
1391c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
13924b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
13934b27c08aSLois Curfman McInnes -    f - 2-norm of function
13949b94acceSBarry Smith 
139536851e7fSLois Curfman McInnes    Level: advanced
139636851e7fSLois Curfman McInnes 
13979b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
13989b94acceSBarry Smith 
13994b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR()
14009b94acceSBarry Smith @*/
140163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
14029b94acceSBarry Smith {
14033a40ed3dSBarry Smith   PetscFunctionBegin;
14044482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14059b94acceSBarry Smith   (snes)->converged = func;
14069b94acceSBarry Smith   (snes)->cnvP      = cctx;
14073a40ed3dSBarry Smith   PetscFunctionReturn(0);
14089b94acceSBarry Smith }
14099b94acceSBarry Smith 
14104a2ae208SSatish Balay #undef __FUNCT__
14114a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
141252baeb72SSatish Balay /*@
1413184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1414184914b5SBarry Smith 
1415184914b5SBarry Smith    Not Collective
1416184914b5SBarry Smith 
1417184914b5SBarry Smith    Input Parameter:
1418184914b5SBarry Smith .  snes - the SNES context
1419184914b5SBarry Smith 
1420184914b5SBarry Smith    Output Parameter:
1421e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1422184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1423184914b5SBarry Smith 
1424184914b5SBarry Smith    Level: intermediate
1425184914b5SBarry Smith 
1426184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1427184914b5SBarry Smith 
1428184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1429184914b5SBarry Smith 
14304b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1431184914b5SBarry Smith @*/
143263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1433184914b5SBarry Smith {
1434184914b5SBarry Smith   PetscFunctionBegin;
14354482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14364482741eSBarry Smith   PetscValidPointer(reason,2);
1437184914b5SBarry Smith   *reason = snes->reason;
1438184914b5SBarry Smith   PetscFunctionReturn(0);
1439184914b5SBarry Smith }
1440184914b5SBarry Smith 
14414a2ae208SSatish Balay #undef __FUNCT__
14424a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1443c9005455SLois Curfman McInnes /*@
1444c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1445c9005455SLois Curfman McInnes 
1446fee21e36SBarry Smith    Collective on SNES
1447fee21e36SBarry Smith 
1448c7afd0dbSLois Curfman McInnes    Input Parameters:
1449c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1450c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1451cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1452758f92a0SBarry Smith .  na  - size of a and its
145364731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1454758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1455c7afd0dbSLois Curfman McInnes 
1456c9005455SLois Curfman McInnes    Notes:
14574b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1458c9005455SLois Curfman McInnes    at each step.
1459c9005455SLois Curfman McInnes 
1460c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1461c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1462c9005455SLois Curfman McInnes    during the section of code that is being timed.
1463c9005455SLois Curfman McInnes 
146436851e7fSLois Curfman McInnes    Level: intermediate
146536851e7fSLois Curfman McInnes 
1466c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1467758f92a0SBarry Smith 
146808405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1469758f92a0SBarry Smith 
1470c9005455SLois Curfman McInnes @*/
147163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset)
1472c9005455SLois Curfman McInnes {
14733a40ed3dSBarry Smith   PetscFunctionBegin;
14744482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14754482741eSBarry Smith   if (na) PetscValidScalarPointer(a,2);
1476c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1477758f92a0SBarry Smith   snes->conv_hist_its   = its;
1478758f92a0SBarry Smith   snes->conv_hist_max   = na;
1479758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1480758f92a0SBarry Smith   PetscFunctionReturn(0);
1481758f92a0SBarry Smith }
1482758f92a0SBarry Smith 
14834a2ae208SSatish Balay #undef __FUNCT__
14844a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
14850c4c9dddSBarry Smith /*@C
1486758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1487758f92a0SBarry Smith 
1488758f92a0SBarry Smith    Collective on SNES
1489758f92a0SBarry Smith 
1490758f92a0SBarry Smith    Input Parameter:
1491758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1492758f92a0SBarry Smith 
1493758f92a0SBarry Smith    Output Parameters:
1494758f92a0SBarry Smith .  a   - array to hold history
1495758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1496758f92a0SBarry Smith          negative if not converged) for each solve.
1497758f92a0SBarry Smith -  na  - size of a and its
1498758f92a0SBarry Smith 
1499758f92a0SBarry Smith    Notes:
1500758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1501758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1502758f92a0SBarry Smith 
1503758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1504758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1505758f92a0SBarry Smith    during the section of code that is being timed.
1506758f92a0SBarry Smith 
1507758f92a0SBarry Smith    Level: intermediate
1508758f92a0SBarry Smith 
1509758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1510758f92a0SBarry Smith 
1511758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1512758f92a0SBarry Smith 
1513758f92a0SBarry Smith @*/
151463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na)
1515758f92a0SBarry Smith {
1516758f92a0SBarry Smith   PetscFunctionBegin;
15174482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1518758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1519758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1520758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
15213a40ed3dSBarry Smith   PetscFunctionReturn(0);
1522c9005455SLois Curfman McInnes }
1523c9005455SLois Curfman McInnes 
1524e74ef692SMatthew Knepley #undef __FUNCT__
1525e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1526ac226902SBarry Smith /*@C
152776b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
152876b2cf59SMatthew Knepley   at the beginning of every step of the iteration.
152976b2cf59SMatthew Knepley 
153076b2cf59SMatthew Knepley   Collective on SNES
153176b2cf59SMatthew Knepley 
153276b2cf59SMatthew Knepley   Input Parameters:
153376b2cf59SMatthew Knepley . snes - The nonlinear solver context
153476b2cf59SMatthew Knepley . func - The function
153576b2cf59SMatthew Knepley 
153676b2cf59SMatthew Knepley   Calling sequence of func:
1537b5d30489SBarry Smith . func (SNES snes, PetscInt step);
153876b2cf59SMatthew Knepley 
153976b2cf59SMatthew Knepley . step - The current step of the iteration
154076b2cf59SMatthew Knepley 
154176b2cf59SMatthew Knepley   Level: intermediate
154276b2cf59SMatthew Knepley 
154376b2cf59SMatthew Knepley .keywords: SNES, update
1544b5d30489SBarry Smith 
154576b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
154676b2cf59SMatthew Knepley @*/
154763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt))
154876b2cf59SMatthew Knepley {
154976b2cf59SMatthew Knepley   PetscFunctionBegin;
15504482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
155176b2cf59SMatthew Knepley   snes->update = func;
155276b2cf59SMatthew Knepley   PetscFunctionReturn(0);
155376b2cf59SMatthew Knepley }
155476b2cf59SMatthew Knepley 
1555e74ef692SMatthew Knepley #undef __FUNCT__
1556e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
155776b2cf59SMatthew Knepley /*@
155876b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
155976b2cf59SMatthew Knepley 
156076b2cf59SMatthew Knepley   Not collective
156176b2cf59SMatthew Knepley 
156276b2cf59SMatthew Knepley   Input Parameters:
156376b2cf59SMatthew Knepley . snes - The nonlinear solver context
156476b2cf59SMatthew Knepley . step - The current step of the iteration
156576b2cf59SMatthew Knepley 
1566205452f4SMatthew Knepley   Level: intermediate
1567205452f4SMatthew Knepley 
156876b2cf59SMatthew Knepley .keywords: SNES, update
156976b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
157076b2cf59SMatthew Knepley @*/
157163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step)
157276b2cf59SMatthew Knepley {
157376b2cf59SMatthew Knepley   PetscFunctionBegin;
157476b2cf59SMatthew Knepley   PetscFunctionReturn(0);
157576b2cf59SMatthew Knepley }
157676b2cf59SMatthew Knepley 
15774a2ae208SSatish Balay #undef __FUNCT__
15784a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
15799b94acceSBarry Smith /*
15809b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
15819b94acceSBarry Smith    positive parameter delta.
15829b94acceSBarry Smith 
15839b94acceSBarry Smith     Input Parameters:
1584c7afd0dbSLois Curfman McInnes +   snes - the SNES context
15859b94acceSBarry Smith .   y - approximate solution of linear system
15869b94acceSBarry Smith .   fnorm - 2-norm of current function
1587c7afd0dbSLois Curfman McInnes -   delta - trust region size
15889b94acceSBarry Smith 
15899b94acceSBarry Smith     Output Parameters:
1590c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
15919b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
15929b94acceSBarry Smith     region, and exceeds zero otherwise.
1593c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
15949b94acceSBarry Smith 
15959b94acceSBarry Smith     Note:
15964b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
15979b94acceSBarry Smith     is set to be the maximum allowable step size.
15989b94acceSBarry Smith 
15999b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
16009b94acceSBarry Smith */
1601dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
16029b94acceSBarry Smith {
1603064f8208SBarry Smith   PetscReal      nrm;
1604ea709b57SSatish Balay   PetscScalar    cnorm;
1605dfbe8321SBarry Smith   PetscErrorCode ierr;
16063a40ed3dSBarry Smith 
16073a40ed3dSBarry Smith   PetscFunctionBegin;
16084482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16094482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1610c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1611184914b5SBarry Smith 
1612064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1613064f8208SBarry Smith   if (nrm > *delta) {
1614064f8208SBarry Smith      nrm = *delta/nrm;
1615064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1616064f8208SBarry Smith      cnorm = nrm;
16172dcb1b2aSMatthew Knepley      ierr = VecScale(y,cnorm);CHKERRQ(ierr);
16189b94acceSBarry Smith      *ynorm = *delta;
16199b94acceSBarry Smith   } else {
16209b94acceSBarry Smith      *gpnorm = 0.0;
1621064f8208SBarry Smith      *ynorm = nrm;
16229b94acceSBarry Smith   }
16233a40ed3dSBarry Smith   PetscFunctionReturn(0);
16249b94acceSBarry Smith }
16259b94acceSBarry Smith 
16264a2ae208SSatish Balay #undef __FUNCT__
16274a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
16289b94acceSBarry Smith /*@
1629f69a0ea3SMatthew Knepley    SNESSolve - Solves a nonlinear system F(x) = b.
1630f69a0ea3SMatthew Knepley    Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX().
16319b94acceSBarry Smith 
1632c7afd0dbSLois Curfman McInnes    Collective on SNES
1633c7afd0dbSLois Curfman McInnes 
1634b2002411SLois Curfman McInnes    Input Parameters:
1635c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1636f69a0ea3SMatthew Knepley .  b - the constant part of the equation, or PETSC_NULL to use zero.
163770e92668SMatthew Knepley -  x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution()
16389b94acceSBarry Smith 
1639b2002411SLois Curfman McInnes    Notes:
16408ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
16418ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
16428ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
16438ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
16448ddd3da0SLois Curfman McInnes 
164536851e7fSLois Curfman McInnes    Level: beginner
164636851e7fSLois Curfman McInnes 
16479b94acceSBarry Smith .keywords: SNES, nonlinear, solve
16489b94acceSBarry Smith 
164970e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution()
16509b94acceSBarry Smith @*/
1651f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x)
16529b94acceSBarry Smith {
1653dfbe8321SBarry Smith   PetscErrorCode ierr;
1654f1af5d2fSBarry Smith   PetscTruth     flg;
1655052efed2SBarry Smith 
16563a40ed3dSBarry Smith   PetscFunctionBegin;
16574482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16581302d50aSBarry Smith   if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
165974637425SBarry Smith 
1660f69a0ea3SMatthew Knepley   if (b) {
1661f69a0ea3SMatthew Knepley     ierr = SNESSetRhs(snes, b); CHKERRQ(ierr);
16621096aae1SMatthew Knepley     if (!snes->vec_func) {
16631096aae1SMatthew Knepley       Vec r;
16641096aae1SMatthew Knepley 
16651096aae1SMatthew Knepley       ierr = VecDuplicate(b, &r); CHKERRQ(ierr);
16661096aae1SMatthew Knepley       ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr);
16671096aae1SMatthew Knepley     }
1668f69a0ea3SMatthew Knepley   }
166970e92668SMatthew Knepley   if (x) {
1670f69a0ea3SMatthew Knepley     PetscValidHeaderSpecific(x,VEC_COOKIE,3);
1671f69a0ea3SMatthew Knepley     PetscCheckSameComm(snes,1,x,3);
167270e92668SMatthew Knepley   } else {
167370e92668SMatthew Knepley     ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr);
167470e92668SMatthew Knepley     if (!x) {
167570e92668SMatthew Knepley       ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr);
167670e92668SMatthew Knepley     }
167770e92668SMatthew Knepley   }
167870e92668SMatthew Knepley   snes->vec_sol = snes->vec_sol_always = x;
167970e92668SMatthew Knepley   if (!snes->setupcalled) {
168070e92668SMatthew Knepley     ierr = SNESSetUp(snes);CHKERRQ(ierr);
168170e92668SMatthew Knepley   }
1682abc0a331SBarry Smith   if (snes->conv_hist_reset) snes->conv_hist_len = 0;
1683d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
168450ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1685d5e45103SBarry Smith 
1686d5e45103SBarry Smith   ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN);
1687d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
168838f152feSBarry Smith     /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */
168919717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr);
1690d5e45103SBarry Smith   } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) {
1691d5e45103SBarry Smith     /* translate exception into SNES not converged reason */
1692d5e45103SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
169338f152feSBarry Smith     ierr = 0;
169438f152feSBarry Smith   }
169538f152feSBarry Smith   CHKERRQ(ierr);
1696d5e45103SBarry Smith 
1697d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1698b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
16998b179ff8SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); }
1700da9b6338SBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1701da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
17025968eb51SBarry Smith   if (snes->printreason) {
17035968eb51SBarry Smith     if (snes->reason > 0) {
17049dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17055968eb51SBarry Smith     } else {
17069dcbbd2bSBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr);
17075968eb51SBarry Smith     }
17085968eb51SBarry Smith   }
17095968eb51SBarry Smith 
17103a40ed3dSBarry Smith   PetscFunctionReturn(0);
17119b94acceSBarry Smith }
17129b94acceSBarry Smith 
17139b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
17149b94acceSBarry Smith 
17154a2ae208SSatish Balay #undef __FUNCT__
17164a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
171782bf6240SBarry Smith /*@C
17184b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
17199b94acceSBarry Smith 
1720fee21e36SBarry Smith    Collective on SNES
1721fee21e36SBarry Smith 
1722c7afd0dbSLois Curfman McInnes    Input Parameters:
1723c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1724454a90a3SBarry Smith -  type - a known method
1725c7afd0dbSLois Curfman McInnes 
1726c7afd0dbSLois Curfman McInnes    Options Database Key:
1727454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1728c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1729ae12b187SLois Curfman McInnes 
17309b94acceSBarry Smith    Notes:
1731e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
17324b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
1733c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17344b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
1735c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
17369b94acceSBarry Smith 
1737ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1738ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1739ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1740ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1741ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1742ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1743ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1744ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1745ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1746b0a32e0cSBarry Smith   appropriate method.
174736851e7fSLois Curfman McInnes 
174836851e7fSLois Curfman McInnes   Level: intermediate
1749a703fe33SLois Curfman McInnes 
1750454a90a3SBarry Smith .keywords: SNES, set, type
1751435da068SBarry Smith 
1752435da068SBarry Smith .seealso: SNESType, SNESCreate()
1753435da068SBarry Smith 
17549b94acceSBarry Smith @*/
1755e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type)
17569b94acceSBarry Smith {
1757dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
17586831982aSBarry Smith   PetscTruth     match;
17593a40ed3dSBarry Smith 
17603a40ed3dSBarry Smith   PetscFunctionBegin;
17614482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17624482741eSBarry Smith   PetscValidCharPointer(type,2);
176382bf6240SBarry Smith 
17646831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
17650f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
176682bf6240SBarry Smith 
176782bf6240SBarry Smith   if (snes->setupcalled) {
17684dc4c822SBarry Smith     snes->setupcalled = PETSC_FALSE;
1769e1311b90SBarry Smith     ierr              = (*(snes)->destroy)(snes);CHKERRQ(ierr);
177082bf6240SBarry Smith     snes->data        = 0;
177182bf6240SBarry Smith   }
17727d1a2b2bSBarry Smith 
17739b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
177482bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1775b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1776958c9bccSBarry Smith   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
1777606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
177882bf6240SBarry Smith   snes->data = 0;
17793a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1780454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
17813a40ed3dSBarry Smith   PetscFunctionReturn(0);
17829b94acceSBarry Smith }
17839b94acceSBarry Smith 
1784a847f771SSatish Balay 
17859b94acceSBarry Smith /* --------------------------------------------------------------------- */
17864a2ae208SSatish Balay #undef __FUNCT__
17874a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
178852baeb72SSatish Balay /*@
17899b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1790f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
17919b94acceSBarry Smith 
1792fee21e36SBarry Smith    Not Collective
1793fee21e36SBarry Smith 
179436851e7fSLois Curfman McInnes    Level: advanced
179536851e7fSLois Curfman McInnes 
17969b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
17979b94acceSBarry Smith 
17989b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
17999b94acceSBarry Smith @*/
180063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void)
18019b94acceSBarry Smith {
1802dfbe8321SBarry Smith   PetscErrorCode ierr;
180382bf6240SBarry Smith 
18043a40ed3dSBarry Smith   PetscFunctionBegin;
180582bf6240SBarry Smith   if (SNESList) {
1806b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
180782bf6240SBarry Smith     SNESList = 0;
18089b94acceSBarry Smith   }
18094c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
18103a40ed3dSBarry Smith   PetscFunctionReturn(0);
18119b94acceSBarry Smith }
18129b94acceSBarry Smith 
18134a2ae208SSatish Balay #undef __FUNCT__
18144a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
18159b94acceSBarry Smith /*@C
18169a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
18179b94acceSBarry Smith 
1818c7afd0dbSLois Curfman McInnes    Not Collective
1819c7afd0dbSLois Curfman McInnes 
18209b94acceSBarry Smith    Input Parameter:
18214b0e389bSBarry Smith .  snes - nonlinear solver context
18229b94acceSBarry Smith 
18239b94acceSBarry Smith    Output Parameter:
18243a7fca6bSBarry Smith .  type - SNES method (a character string)
18259b94acceSBarry Smith 
182636851e7fSLois Curfman McInnes    Level: intermediate
182736851e7fSLois Curfman McInnes 
1828454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
18299b94acceSBarry Smith @*/
183063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type)
18319b94acceSBarry Smith {
18323a40ed3dSBarry Smith   PetscFunctionBegin;
18334482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18344482741eSBarry Smith   PetscValidPointer(type,2);
1835454a90a3SBarry Smith   *type = snes->type_name;
18363a40ed3dSBarry Smith   PetscFunctionReturn(0);
18379b94acceSBarry Smith }
18389b94acceSBarry Smith 
18394a2ae208SSatish Balay #undef __FUNCT__
18404a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
184152baeb72SSatish Balay /*@
18429b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
18439b94acceSBarry Smith    stored.
18449b94acceSBarry Smith 
1845c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1846c7afd0dbSLois Curfman McInnes 
18479b94acceSBarry Smith    Input Parameter:
18489b94acceSBarry Smith .  snes - the SNES context
18499b94acceSBarry Smith 
18509b94acceSBarry Smith    Output Parameter:
18519b94acceSBarry Smith .  x - the solution
18529b94acceSBarry Smith 
185370e92668SMatthew Knepley    Level: intermediate
185436851e7fSLois Curfman McInnes 
18559b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
18569b94acceSBarry Smith 
185770e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
18589b94acceSBarry Smith @*/
185963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x)
18609b94acceSBarry Smith {
18613a40ed3dSBarry Smith   PetscFunctionBegin;
18624482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18634482741eSBarry Smith   PetscValidPointer(x,2);
18649b94acceSBarry Smith   *x = snes->vec_sol_always;
18653a40ed3dSBarry Smith   PetscFunctionReturn(0);
18669b94acceSBarry Smith }
18679b94acceSBarry Smith 
18684a2ae208SSatish Balay #undef __FUNCT__
186970e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution"
187070e92668SMatthew Knepley /*@C
187170e92668SMatthew Knepley    SNESSetSolution - Sets the vector where the approximate solution is stored.
187270e92668SMatthew Knepley 
187370e92668SMatthew Knepley    Not Collective, but Vec is parallel if SNES is parallel
187470e92668SMatthew Knepley 
187570e92668SMatthew Knepley    Input Parameters:
187670e92668SMatthew Knepley +  snes - the SNES context
187770e92668SMatthew Knepley -  x - the solution
187870e92668SMatthew Knepley 
187970e92668SMatthew Knepley    Output Parameter:
188070e92668SMatthew Knepley 
188170e92668SMatthew Knepley    Level: intermediate
188270e92668SMatthew Knepley 
188342219521SBarry Smith    Notes: this is not normally used, rather one simply calls SNESSolve() with
188442219521SBarry Smith           the appropriate solution vector.
188542219521SBarry Smith 
188670e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution
188770e92668SMatthew Knepley 
188870e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate()
188970e92668SMatthew Knepley @*/
189070e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x)
189170e92668SMatthew Knepley {
189270e92668SMatthew Knepley   PetscFunctionBegin;
189370e92668SMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
189470e92668SMatthew Knepley   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
189570e92668SMatthew Knepley   PetscCheckSameComm(snes,1,x,2);
189670e92668SMatthew Knepley   snes->vec_sol_always = x;
189770e92668SMatthew Knepley   PetscFunctionReturn(0);
189870e92668SMatthew Knepley }
189970e92668SMatthew Knepley 
190070e92668SMatthew Knepley #undef __FUNCT__
19014a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
190252baeb72SSatish Balay /*@
19039b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
19049b94acceSBarry Smith    stored.
19059b94acceSBarry Smith 
1906c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1907c7afd0dbSLois Curfman McInnes 
19089b94acceSBarry Smith    Input Parameter:
19099b94acceSBarry Smith .  snes - the SNES context
19109b94acceSBarry Smith 
19119b94acceSBarry Smith    Output Parameter:
19129b94acceSBarry Smith .  x - the solution update
19139b94acceSBarry Smith 
191436851e7fSLois Curfman McInnes    Level: advanced
191536851e7fSLois Curfman McInnes 
19169b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
19179b94acceSBarry Smith 
19189b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
19199b94acceSBarry Smith @*/
192063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x)
19219b94acceSBarry Smith {
19223a40ed3dSBarry Smith   PetscFunctionBegin;
19234482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19244482741eSBarry Smith   PetscValidPointer(x,2);
19259b94acceSBarry Smith   *x = snes->vec_sol_update_always;
19263a40ed3dSBarry Smith   PetscFunctionReturn(0);
19279b94acceSBarry Smith }
19289b94acceSBarry Smith 
19294a2ae208SSatish Balay #undef __FUNCT__
19304a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
19319b94acceSBarry Smith /*@C
19323638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19339b94acceSBarry Smith 
1934c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1935c7afd0dbSLois Curfman McInnes 
19369b94acceSBarry Smith    Input Parameter:
19379b94acceSBarry Smith .  snes - the SNES context
19389b94acceSBarry Smith 
19399b94acceSBarry Smith    Output Parameter:
19407bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
194170e92668SMatthew Knepley .  func - the function (or PETSC_NULL)
194270e92668SMatthew Knepley -  ctx - the function context (or PETSC_NULL)
19439b94acceSBarry Smith 
194436851e7fSLois Curfman McInnes    Level: advanced
194536851e7fSLois Curfman McInnes 
1946a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19479b94acceSBarry Smith 
19484b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
19499b94acceSBarry Smith @*/
195070e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx)
19519b94acceSBarry Smith {
19523a40ed3dSBarry Smith   PetscFunctionBegin;
19534482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19547bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
195500036973SBarry Smith   if (func) *func = snes->computefunction;
195670e92668SMatthew Knepley   if (ctx)  *ctx  = snes->funP;
19573a40ed3dSBarry Smith   PetscFunctionReturn(0);
19589b94acceSBarry Smith }
19599b94acceSBarry Smith 
19604a2ae208SSatish Balay #undef __FUNCT__
19614a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
19623c7409f5SSatish Balay /*@C
19633c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
1964d850072dSLois Curfman McInnes    SNES options in the database.
19653c7409f5SSatish Balay 
1966fee21e36SBarry Smith    Collective on SNES
1967fee21e36SBarry Smith 
1968c7afd0dbSLois Curfman McInnes    Input Parameter:
1969c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1970c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
1971c7afd0dbSLois Curfman McInnes 
1972d850072dSLois Curfman McInnes    Notes:
1973a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1974c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
1975d850072dSLois Curfman McInnes 
197636851e7fSLois Curfman McInnes    Level: advanced
197736851e7fSLois Curfman McInnes 
19783c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
1979a86d99e1SLois Curfman McInnes 
1980a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
19813c7409f5SSatish Balay @*/
198263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[])
19833c7409f5SSatish Balay {
1984dfbe8321SBarry Smith   PetscErrorCode ierr;
19853c7409f5SSatish Balay 
19863a40ed3dSBarry Smith   PetscFunctionBegin;
19874482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1988639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
198994b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
19903a40ed3dSBarry Smith   PetscFunctionReturn(0);
19913c7409f5SSatish Balay }
19923c7409f5SSatish Balay 
19934a2ae208SSatish Balay #undef __FUNCT__
19944a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
19953c7409f5SSatish Balay /*@C
1996f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
1997d850072dSLois Curfman McInnes    SNES options in the database.
19983c7409f5SSatish Balay 
1999fee21e36SBarry Smith    Collective on SNES
2000fee21e36SBarry Smith 
2001c7afd0dbSLois Curfman McInnes    Input Parameters:
2002c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2003c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2004c7afd0dbSLois Curfman McInnes 
2005d850072dSLois Curfman McInnes    Notes:
2006a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2007c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2008d850072dSLois Curfman McInnes 
200936851e7fSLois Curfman McInnes    Level: advanced
201036851e7fSLois Curfman McInnes 
20113c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2012a86d99e1SLois Curfman McInnes 
2013a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20143c7409f5SSatish Balay @*/
201563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[])
20163c7409f5SSatish Balay {
2017dfbe8321SBarry Smith   PetscErrorCode ierr;
20183c7409f5SSatish Balay 
20193a40ed3dSBarry Smith   PetscFunctionBegin;
20204482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2021639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
202294b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20233a40ed3dSBarry Smith   PetscFunctionReturn(0);
20243c7409f5SSatish Balay }
20253c7409f5SSatish Balay 
20264a2ae208SSatish Balay #undef __FUNCT__
20274a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
20289ab63eb5SSatish Balay /*@C
20293c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20303c7409f5SSatish Balay    SNES options in the database.
20313c7409f5SSatish Balay 
2032c7afd0dbSLois Curfman McInnes    Not Collective
2033c7afd0dbSLois Curfman McInnes 
20343c7409f5SSatish Balay    Input Parameter:
20353c7409f5SSatish Balay .  snes - the SNES context
20363c7409f5SSatish Balay 
20373c7409f5SSatish Balay    Output Parameter:
20383c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20393c7409f5SSatish Balay 
20409ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
20419ab63eb5SSatish Balay    sufficient length to hold the prefix.
20429ab63eb5SSatish Balay 
204336851e7fSLois Curfman McInnes    Level: advanced
204436851e7fSLois Curfman McInnes 
20453c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2046a86d99e1SLois Curfman McInnes 
2047a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20483c7409f5SSatish Balay @*/
2049e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[])
20503c7409f5SSatish Balay {
2051dfbe8321SBarry Smith   PetscErrorCode ierr;
20523c7409f5SSatish Balay 
20533a40ed3dSBarry Smith   PetscFunctionBegin;
20544482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2055639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
20563a40ed3dSBarry Smith   PetscFunctionReturn(0);
20573c7409f5SSatish Balay }
20583c7409f5SSatish Balay 
2059b2002411SLois Curfman McInnes 
20604a2ae208SSatish Balay #undef __FUNCT__
20614a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
20623cea93caSBarry Smith /*@C
20633cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
20643cea93caSBarry Smith 
20657f6c08e0SMatthew Knepley   Level: advanced
20663cea93caSBarry Smith @*/
206763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2068b2002411SLois Curfman McInnes {
2069e2d1d2b7SBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
2070dfbe8321SBarry Smith   PetscErrorCode ierr;
2071b2002411SLois Curfman McInnes 
2072b2002411SLois Curfman McInnes   PetscFunctionBegin;
2073b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2074c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2075b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2076b2002411SLois Curfman McInnes }
2077da9b6338SBarry Smith 
2078da9b6338SBarry Smith #undef __FUNCT__
2079da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
208063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes)
2081da9b6338SBarry Smith {
2082dfbe8321SBarry Smith   PetscErrorCode ierr;
208377431f27SBarry Smith   PetscInt       N,i,j;
2084da9b6338SBarry Smith   Vec            u,uh,fh;
2085da9b6338SBarry Smith   PetscScalar    value;
2086da9b6338SBarry Smith   PetscReal      norm;
2087da9b6338SBarry Smith 
2088da9b6338SBarry Smith   PetscFunctionBegin;
2089da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2090da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2091da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2092da9b6338SBarry Smith 
2093da9b6338SBarry Smith   /* currently only works for sequential */
2094da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2095da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2096da9b6338SBarry Smith   for (i=0; i<N; i++) {
2097da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
209877431f27SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr);
2099da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2100ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2101da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
21023ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2103da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
210477431f27SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %D %18.16e\n",j,norm);CHKERRQ(ierr);
2105da9b6338SBarry Smith       value = -value;
2106da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2107da9b6338SBarry Smith     }
2108da9b6338SBarry Smith   }
2109da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2110da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2111da9b6338SBarry Smith   PetscFunctionReturn(0);
2112da9b6338SBarry Smith }
2113