xref: /petsc/src/snes/interface/snes.c (revision 32a4b47a72620a1b77aa3c4e8ddc8668bd2d097b)
173f4d377SMatthew Knepley /*$Id: snes.c,v 1.235 2001/08/21 21:03:49 bsmith Exp $*/
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 */
98ba1e511SMatthew Knepley int SNES_COOKIE;
108ba1e511SMatthew Knepley int MATSNESMFCTX_COOKIE;
11d5ba7fb7SMatthew Knepley int SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval, SNES_MinimizationFunctionEval, SNES_GradientEval;
12d5ba7fb7SMatthew Knepley int SNES_HessianEval;
1382bf6240SBarry Smith 
144a2ae208SSatish Balay #undef __FUNCT__
15a09944afSBarry Smith #define __FUNCT__ "SNESGetProblemType"
16a09944afSBarry Smith /*@C
17a09944afSBarry Smith    SNESGetProblemType -Indicates if SNES is solving a nonlinear system or a minimization
18a09944afSBarry Smith 
19a09944afSBarry Smith    Not Collective
20a09944afSBarry Smith 
21a09944afSBarry Smith    Input Parameter:
22a09944afSBarry Smith .  SNES - the SNES context
23a09944afSBarry Smith 
24a09944afSBarry Smith    Output Parameter:
25a09944afSBarry Smith .   type - SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
26a09944afSBarry Smith    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
27a09944afSBarry Smith 
28a09944afSBarry Smith    Level: intermediate
29a09944afSBarry Smith 
30a09944afSBarry Smith .keywords: SNES, problem type
31a09944afSBarry Smith 
32a09944afSBarry Smith .seealso: SNESCreate()
33a09944afSBarry Smith @*/
34a09944afSBarry Smith int SNESGetProblemType(SNES snes,SNESProblemType *type)
35a09944afSBarry Smith {
36a09944afSBarry Smith   PetscFunctionBegin;
37a09944afSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
38a09944afSBarry Smith   *type = snes->method_class;
39a09944afSBarry Smith   PetscFunctionReturn(0);
40a09944afSBarry Smith }
41a09944afSBarry Smith 
42a09944afSBarry Smith #undef __FUNCT__
434a2ae208SSatish Balay #define __FUNCT__ "SNESView"
447e2c5f70SBarry Smith /*@C
459b94acceSBarry Smith    SNESView - Prints the SNES data structure.
469b94acceSBarry Smith 
474c49b128SBarry Smith    Collective on SNES
48fee21e36SBarry Smith 
49c7afd0dbSLois Curfman McInnes    Input Parameters:
50c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
51c7afd0dbSLois Curfman McInnes -  viewer - visualization context
52c7afd0dbSLois Curfman McInnes 
539b94acceSBarry Smith    Options Database Key:
54c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
559b94acceSBarry Smith 
569b94acceSBarry Smith    Notes:
579b94acceSBarry Smith    The available visualization contexts include
58b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
59b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
60c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
61c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
62c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
639b94acceSBarry Smith 
643e081fefSLois Curfman McInnes    The user can open an alternative visualization context with
65b0a32e0cSBarry Smith    PetscViewerASCIIOpen() - output to a specified file.
669b94acceSBarry Smith 
6736851e7fSLois Curfman McInnes    Level: beginner
6836851e7fSLois Curfman McInnes 
699b94acceSBarry Smith .keywords: SNES, view
709b94acceSBarry Smith 
71b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
729b94acceSBarry Smith @*/
73b0a32e0cSBarry Smith int SNESView(SNES snes,PetscViewer viewer)
749b94acceSBarry Smith {
759b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
769b94acceSBarry Smith   int                 ierr;
779b94acceSBarry Smith   SLES                sles;
78454a90a3SBarry Smith   char                *type;
796831982aSBarry Smith   PetscTruth          isascii,isstring;
809b94acceSBarry Smith 
813a40ed3dSBarry Smith   PetscFunctionBegin;
8274679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
83b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);
84b0a32e0cSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE);
856831982aSBarry Smith   PetscCheckSameComm(snes,viewer);
8674679c65SBarry Smith 
87b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
88b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
890f5bd95cSBarry Smith   if (isascii) {
903a7fca6bSBarry Smith     if (snes->prefix) {
913a7fca6bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr);
923a7fca6bSBarry Smith     } else {
93b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
943a7fca6bSBarry Smith     }
95454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
96454a90a3SBarry Smith     if (type) {
97b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
98184914b5SBarry Smith     } else {
99b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
100184914b5SBarry Smith     }
1010ef38995SBarry Smith     if (snes->view) {
102b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1030ef38995SBarry Smith       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
104b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1050ef38995SBarry Smith     }
106b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
10777d8c4bbSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",
10877d8c4bbSBarry Smith                  snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr);
109b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
110b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
1110ef38995SBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
112b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr);
1130ef38995SBarry Smith     }
1149b94acceSBarry Smith     if (snes->ksp_ewconv) {
1159b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
1169b94acceSBarry Smith       if (kctx) {
117b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
118b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
119b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
1209b94acceSBarry Smith       }
1219b94acceSBarry Smith     }
1220f5bd95cSBarry Smith   } else if (isstring) {
123454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
124b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
12519bcc07fSBarry Smith   }
12677ed5343SBarry Smith   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
127b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1289b94acceSBarry Smith   ierr = SLESView(sles,viewer);CHKERRQ(ierr);
129b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1319b94acceSBarry Smith }
1329b94acceSBarry Smith 
13376b2cf59SMatthew Knepley /*
13476b2cf59SMatthew Knepley   We retain a list of functions that also take SNES command
13576b2cf59SMatthew Knepley   line options. These are called at the end SNESSetFromOptions()
13676b2cf59SMatthew Knepley */
13776b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5
13876b2cf59SMatthew Knepley static int numberofsetfromoptions;
13976b2cf59SMatthew Knepley static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
14076b2cf59SMatthew Knepley 
141e74ef692SMatthew Knepley #undef __FUNCT__
142e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker"
14376b2cf59SMatthew Knepley /*@
14476b2cf59SMatthew Knepley   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
14576b2cf59SMatthew Knepley 
14676b2cf59SMatthew Knepley   Not Collective
14776b2cf59SMatthew Knepley 
14876b2cf59SMatthew Knepley   Input Parameter:
14976b2cf59SMatthew Knepley . snescheck - function that checks for options
15076b2cf59SMatthew Knepley 
15176b2cf59SMatthew Knepley   Level: developer
15276b2cf59SMatthew Knepley 
15376b2cf59SMatthew Knepley .seealso: SNESSetFromOptions()
15476b2cf59SMatthew Knepley @*/
15576b2cf59SMatthew Knepley int SNESAddOptionsChecker(int (*snescheck)(SNES))
15676b2cf59SMatthew Knepley {
15776b2cf59SMatthew Knepley   PetscFunctionBegin;
15876b2cf59SMatthew Knepley   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
15976b2cf59SMatthew Knepley     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS);
16076b2cf59SMatthew Knepley   }
16176b2cf59SMatthew Knepley 
16276b2cf59SMatthew Knepley   othersetfromoptions[numberofsetfromoptions++] = snescheck;
16376b2cf59SMatthew Knepley   PetscFunctionReturn(0);
16476b2cf59SMatthew Knepley }
16576b2cf59SMatthew Knepley 
1664a2ae208SSatish Balay #undef __FUNCT__
1674a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions"
1689b94acceSBarry Smith /*@
169682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
1709b94acceSBarry Smith 
171c7afd0dbSLois Curfman McInnes    Collective on SNES
172c7afd0dbSLois Curfman McInnes 
1739b94acceSBarry Smith    Input Parameter:
1749b94acceSBarry Smith .  snes - the SNES context
1759b94acceSBarry Smith 
17636851e7fSLois Curfman McInnes    Options Database Keys:
1776831982aSBarry Smith +  -snes_type <type> - ls, tr, umls, umtr, test
17882738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
17982738288SBarry Smith                 of the change in the solution between steps
180b39c3a46SLois Curfman McInnes .  -snes_atol <atol> - absolute tolerance of residual norm
181b39c3a46SLois Curfman McInnes .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
182b39c3a46SLois Curfman McInnes .  -snes_max_it <max_it> - maximum number of iterations
183b39c3a46SLois Curfman McInnes .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
18450ffb88aSMatthew Knepley .  -snes_max_fail <max_fail> - maximum number of failures
185b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
18682738288SBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
18782738288SBarry Smith                                solver; hence iterations will continue until max_it
1881fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
18982738288SBarry Smith                                of convergence test
19082738288SBarry Smith .  -snes_monitor - prints residual norm at each iteration
191d132466eSBarry Smith .  -snes_vecmonitor - plots solution at each iteration
192d132466eSBarry Smith .  -snes_vecmonitor_update - plots update to solution at each iteration
19382738288SBarry Smith .  -snes_xmonitor - plots residual norm at each iteration
194e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
19536851e7fSLois Curfman McInnes -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
19682738288SBarry Smith 
19782738288SBarry Smith     Options Database for Eisenstat-Walker method:
19882738288SBarry Smith +  -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence
19982738288SBarry Smith .  -snes_ksp_eq_version ver - version of  Eisenstat-Walker method
20036851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
20136851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
20236851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
20336851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
20436851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
20536851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
20682738288SBarry Smith 
20711ca99fdSLois Curfman McInnes    Notes:
20811ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
20911ca99fdSLois Curfman McInnes    the users manual.
21083e2fdc7SBarry Smith 
21136851e7fSLois Curfman McInnes    Level: beginner
21236851e7fSLois Curfman McInnes 
2139b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
2149b94acceSBarry Smith 
21569ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
2169b94acceSBarry Smith @*/
2179b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
2189b94acceSBarry Smith {
2199b94acceSBarry Smith   SLES                sles;
220186905e3SBarry Smith   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
221f1af5d2fSBarry Smith   PetscTruth          flg;
22276b2cf59SMatthew Knepley   int                 ierr, i;
223186905e3SBarry Smith   char                *deft,type[256];
2249b94acceSBarry Smith 
2253a40ed3dSBarry Smith   PetscFunctionBegin;
22677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
227ca161407SBarry Smith 
228b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
229186905e3SBarry Smith     if (snes->type_name) {
230186905e3SBarry Smith       deft = snes->type_name;
231186905e3SBarry Smith     } else {
232186905e3SBarry Smith       if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
233186905e3SBarry Smith         deft = SNESEQLS;
234186905e3SBarry Smith       } else {
235186905e3SBarry Smith         deft = SNESUMTR;
236d64ed03dSBarry Smith       }
237d64ed03dSBarry Smith     }
2384bbc92c1SBarry Smith 
239186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
240b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
241d64ed03dSBarry Smith     if (flg) {
242186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
243186905e3SBarry Smith     } else if (!snes->type_name) {
244186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
245d64ed03dSBarry Smith     }
246909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
24793c39befSBarry Smith 
24887828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
24987828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
250186905e3SBarry Smith 
25187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
252b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
253b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
25450ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
25587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_fmin","Minimization function tolerance","SNESSetMinimizationFunctionTolerance",snes->fmin,&snes->fmin,0);CHKERRQ(ierr);
256186905e3SBarry Smith 
257b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
258186905e3SBarry Smith 
259b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
26087828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
26187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
26287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
26387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
26487828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
26587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
266186905e3SBarry Smith 
267b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
26893c39befSBarry Smith     if (flg) {snes->converged = 0;}
269b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
2705cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
271b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
272b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
2733a7fca6bSBarry Smith     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
2743a7fca6bSBarry Smith     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
275b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
276b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
277b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
278b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
279b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2807c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
2815ed2d596SBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
2825ed2d596SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
283b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
284186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
285e24b481bSBarry Smith 
286b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2873c7409f5SSatish Balay     if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
288186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
289b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
290981c4779SBarry Smith     } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
291186905e3SBarry Smith       ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeHessian,snes->funP);CHKERRQ(ierr);
292b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2939b94acceSBarry Smith     }
294639f9d9dSBarry Smith 
29576b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
29676b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);                                                             CHKERRQ(ierr);
29776b2cf59SMatthew Knepley     }
29876b2cf59SMatthew Knepley 
299186905e3SBarry Smith     if (snes->setfromoptions) {
300186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
301639f9d9dSBarry Smith     }
302639f9d9dSBarry Smith 
303b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
3044bbc92c1SBarry Smith 
3059b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
3069b94acceSBarry Smith   ierr = SLESSetFromOptions(sles);CHKERRQ(ierr);
30793993e2dSLois Curfman McInnes 
3083a40ed3dSBarry Smith   PetscFunctionReturn(0);
3099b94acceSBarry Smith }
3109b94acceSBarry Smith 
311a847f771SSatish Balay 
3124a2ae208SSatish Balay #undef __FUNCT__
3134a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
3149b94acceSBarry Smith /*@
3159b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3169b94acceSBarry Smith    the nonlinear solvers.
3179b94acceSBarry Smith 
318fee21e36SBarry Smith    Collective on SNES
319fee21e36SBarry Smith 
320c7afd0dbSLois Curfman McInnes    Input Parameters:
321c7afd0dbSLois Curfman McInnes +  snes - the SNES context
322c7afd0dbSLois Curfman McInnes -  usrP - optional user context
323c7afd0dbSLois Curfman McInnes 
32436851e7fSLois Curfman McInnes    Level: intermediate
32536851e7fSLois Curfman McInnes 
3269b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3279b94acceSBarry Smith 
3289b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3299b94acceSBarry Smith @*/
3309b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3319b94acceSBarry Smith {
3323a40ed3dSBarry Smith   PetscFunctionBegin;
33377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3349b94acceSBarry Smith   snes->user		= usrP;
3353a40ed3dSBarry Smith   PetscFunctionReturn(0);
3369b94acceSBarry Smith }
33774679c65SBarry Smith 
3384a2ae208SSatish Balay #undef __FUNCT__
3394a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3409b94acceSBarry Smith /*@C
3419b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3429b94acceSBarry Smith    nonlinear solvers.
3439b94acceSBarry Smith 
344c7afd0dbSLois Curfman McInnes    Not Collective
345c7afd0dbSLois Curfman McInnes 
3469b94acceSBarry Smith    Input Parameter:
3479b94acceSBarry Smith .  snes - SNES context
3489b94acceSBarry Smith 
3499b94acceSBarry Smith    Output Parameter:
3509b94acceSBarry Smith .  usrP - user context
3519b94acceSBarry Smith 
35236851e7fSLois Curfman McInnes    Level: intermediate
35336851e7fSLois Curfman McInnes 
3549b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3559b94acceSBarry Smith 
3569b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3579b94acceSBarry Smith @*/
3589b94acceSBarry Smith int SNESGetApplicationContext(SNES snes,void **usrP)
3599b94acceSBarry Smith {
3603a40ed3dSBarry Smith   PetscFunctionBegin;
36177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3629b94acceSBarry Smith   *usrP = snes->user;
3633a40ed3dSBarry Smith   PetscFunctionReturn(0);
3649b94acceSBarry Smith }
36574679c65SBarry Smith 
3664a2ae208SSatish Balay #undef __FUNCT__
3674a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3689b94acceSBarry Smith /*@
369c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
370c8228a4eSBarry Smith    at this time.
3719b94acceSBarry Smith 
372c7afd0dbSLois Curfman McInnes    Not Collective
373c7afd0dbSLois Curfman McInnes 
3749b94acceSBarry Smith    Input Parameter:
3759b94acceSBarry Smith .  snes - SNES context
3769b94acceSBarry Smith 
3779b94acceSBarry Smith    Output Parameter:
3789b94acceSBarry Smith .  iter - iteration number
3799b94acceSBarry Smith 
380c8228a4eSBarry Smith    Notes:
381c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
382c8228a4eSBarry Smith 
383c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
38408405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
38508405cd6SLois Curfman McInnes .vb
38608405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
38708405cd6SLois Curfman McInnes       if (!(it % 2)) {
38808405cd6SLois Curfman McInnes         [compute Jacobian here]
38908405cd6SLois Curfman McInnes       }
39008405cd6SLois Curfman McInnes .ve
391c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
39208405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
393c8228a4eSBarry Smith 
39436851e7fSLois Curfman McInnes    Level: intermediate
39536851e7fSLois Curfman McInnes 
3969b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3979b94acceSBarry Smith @*/
3989b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3999b94acceSBarry Smith {
4003a40ed3dSBarry Smith   PetscFunctionBegin;
40177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
40274679c65SBarry Smith   PetscValidIntPointer(iter);
4039b94acceSBarry Smith   *iter = snes->iter;
4043a40ed3dSBarry Smith   PetscFunctionReturn(0);
4059b94acceSBarry Smith }
40674679c65SBarry Smith 
4074a2ae208SSatish Balay #undef __FUNCT__
4084a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
4099b94acceSBarry Smith /*@
4109b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
4119b94acceSBarry Smith    with SNESSSetFunction().
4129b94acceSBarry Smith 
413c7afd0dbSLois Curfman McInnes    Collective on SNES
414c7afd0dbSLois Curfman McInnes 
4159b94acceSBarry Smith    Input Parameter:
4169b94acceSBarry Smith .  snes - SNES context
4179b94acceSBarry Smith 
4189b94acceSBarry Smith    Output Parameter:
4199b94acceSBarry Smith .  fnorm - 2-norm of function
4209b94acceSBarry Smith 
4219b94acceSBarry Smith    Note:
4229b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
423a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
424a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
4259b94acceSBarry Smith 
42636851e7fSLois Curfman McInnes    Level: intermediate
42736851e7fSLois Curfman McInnes 
4289b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
429a86d99e1SLois Curfman McInnes 
43008405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
4319b94acceSBarry Smith @*/
43287828ca2SBarry Smith int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
4339b94acceSBarry Smith {
4343a40ed3dSBarry Smith   PetscFunctionBegin;
43577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43674679c65SBarry Smith   PetscValidScalarPointer(fnorm);
43774679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
43829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_NONLINEAR_EQUATIONS only");
43974679c65SBarry Smith   }
4409b94acceSBarry Smith   *fnorm = snes->norm;
4413a40ed3dSBarry Smith   PetscFunctionReturn(0);
4429b94acceSBarry Smith }
44374679c65SBarry Smith 
4444a2ae208SSatish Balay #undef __FUNCT__
4454a2ae208SSatish Balay #define __FUNCT__ "SNESGetGradientNorm"
4469b94acceSBarry Smith /*@
4479b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4489b94acceSBarry Smith    with SNESSSetGradient().
4499b94acceSBarry Smith 
450c7afd0dbSLois Curfman McInnes    Collective on SNES
451c7afd0dbSLois Curfman McInnes 
4529b94acceSBarry Smith    Input Parameter:
4539b94acceSBarry Smith .  snes - SNES context
4549b94acceSBarry Smith 
4559b94acceSBarry Smith    Output Parameter:
4569b94acceSBarry Smith .  fnorm - 2-norm of gradient
4579b94acceSBarry Smith 
4589b94acceSBarry Smith    Note:
4599b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
460a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
461a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4629b94acceSBarry Smith 
46336851e7fSLois Curfman McInnes    Level: intermediate
46436851e7fSLois Curfman McInnes 
4659b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
466a86d99e1SLois Curfman McInnes 
467a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4689b94acceSBarry Smith @*/
46987828ca2SBarry Smith int SNESGetGradientNorm(SNES snes,PetscScalar *gnorm)
4709b94acceSBarry Smith {
4713a40ed3dSBarry Smith   PetscFunctionBegin;
47277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
47374679c65SBarry Smith   PetscValidScalarPointer(gnorm);
47474679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
47529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_UNCONSTRAINED_MINIMIZATION only");
47674679c65SBarry Smith   }
4779b94acceSBarry Smith   *gnorm = snes->norm;
4783a40ed3dSBarry Smith   PetscFunctionReturn(0);
4799b94acceSBarry Smith }
48074679c65SBarry Smith 
4814a2ae208SSatish Balay #undef __FUNCT__
4824a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
4839b94acceSBarry Smith /*@
4849b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4859b94acceSBarry Smith    attempted by the nonlinear solver.
4869b94acceSBarry Smith 
487c7afd0dbSLois Curfman McInnes    Not Collective
488c7afd0dbSLois Curfman McInnes 
4899b94acceSBarry Smith    Input Parameter:
4909b94acceSBarry Smith .  snes - SNES context
4919b94acceSBarry Smith 
4929b94acceSBarry Smith    Output Parameter:
4939b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4949b94acceSBarry Smith 
495c96a6f78SLois Curfman McInnes    Notes:
496c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
497c96a6f78SLois Curfman McInnes 
49836851e7fSLois Curfman McInnes    Level: intermediate
49936851e7fSLois Curfman McInnes 
5009b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
5019b94acceSBarry Smith @*/
5029b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
5039b94acceSBarry Smith {
5043a40ed3dSBarry Smith   PetscFunctionBegin;
50577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
50674679c65SBarry Smith   PetscValidIntPointer(nfails);
50750ffb88aSMatthew Knepley   *nfails = snes->numFailures;
50850ffb88aSMatthew Knepley   PetscFunctionReturn(0);
50950ffb88aSMatthew Knepley }
51050ffb88aSMatthew Knepley 
51150ffb88aSMatthew Knepley #undef __FUNCT__
51250ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
51350ffb88aSMatthew Knepley /*@
51450ffb88aSMatthew Knepley    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
51550ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
51650ffb88aSMatthew Knepley 
51750ffb88aSMatthew Knepley    Not Collective
51850ffb88aSMatthew Knepley 
51950ffb88aSMatthew Knepley    Input Parameters:
52050ffb88aSMatthew Knepley +  snes     - SNES context
52150ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
52250ffb88aSMatthew Knepley 
52350ffb88aSMatthew Knepley    Level: intermediate
52450ffb88aSMatthew Knepley 
52550ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
52650ffb88aSMatthew Knepley @*/
52750ffb88aSMatthew Knepley int SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails)
52850ffb88aSMatthew Knepley {
52950ffb88aSMatthew Knepley   PetscFunctionBegin;
53050ffb88aSMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE);
53150ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
53250ffb88aSMatthew Knepley   PetscFunctionReturn(0);
53350ffb88aSMatthew Knepley }
53450ffb88aSMatthew Knepley 
53550ffb88aSMatthew Knepley #undef __FUNCT__
53650ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
53750ffb88aSMatthew Knepley /*@
53850ffb88aSMatthew Knepley    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
53950ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
54050ffb88aSMatthew Knepley 
54150ffb88aSMatthew Knepley    Not Collective
54250ffb88aSMatthew Knepley 
54350ffb88aSMatthew Knepley    Input Parameter:
54450ffb88aSMatthew Knepley .  snes     - SNES context
54550ffb88aSMatthew Knepley 
54650ffb88aSMatthew Knepley    Output Parameter:
54750ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
54850ffb88aSMatthew Knepley 
54950ffb88aSMatthew Knepley    Level: intermediate
55050ffb88aSMatthew Knepley 
55150ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
55250ffb88aSMatthew Knepley @*/
55350ffb88aSMatthew Knepley int SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails)
55450ffb88aSMatthew Knepley {
55550ffb88aSMatthew Knepley   PetscFunctionBegin;
55650ffb88aSMatthew Knepley   PetscValidHeaderSpecific(snes,SNES_COOKIE);
55750ffb88aSMatthew Knepley   PetscValidIntPointer(maxFails);
55850ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
5593a40ed3dSBarry Smith   PetscFunctionReturn(0);
5609b94acceSBarry Smith }
561a847f771SSatish Balay 
5624a2ae208SSatish Balay #undef __FUNCT__
5634a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations"
564c96a6f78SLois Curfman McInnes /*@
565c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
566c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
567c96a6f78SLois Curfman McInnes 
568c7afd0dbSLois Curfman McInnes    Not Collective
569c7afd0dbSLois Curfman McInnes 
570c96a6f78SLois Curfman McInnes    Input Parameter:
571c96a6f78SLois Curfman McInnes .  snes - SNES context
572c96a6f78SLois Curfman McInnes 
573c96a6f78SLois Curfman McInnes    Output Parameter:
574c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
575c96a6f78SLois Curfman McInnes 
576c96a6f78SLois Curfman McInnes    Notes:
577c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
578c96a6f78SLois Curfman McInnes 
57936851e7fSLois Curfman McInnes    Level: intermediate
58036851e7fSLois Curfman McInnes 
581c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
582c96a6f78SLois Curfman McInnes @*/
58386bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
584c96a6f78SLois Curfman McInnes {
5853a40ed3dSBarry Smith   PetscFunctionBegin;
586c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
587c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
588c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5893a40ed3dSBarry Smith   PetscFunctionReturn(0);
590c96a6f78SLois Curfman McInnes }
591c96a6f78SLois Curfman McInnes 
5924a2ae208SSatish Balay #undef __FUNCT__
5934a2ae208SSatish Balay #define __FUNCT__ "SNESGetSLES"
5949b94acceSBarry Smith /*@C
5959b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
5969b94acceSBarry Smith 
597c7afd0dbSLois Curfman McInnes    Not Collective, but if SNES object is parallel, then SLES object is parallel
598c7afd0dbSLois Curfman McInnes 
5999b94acceSBarry Smith    Input Parameter:
6009b94acceSBarry Smith .  snes - the SNES context
6019b94acceSBarry Smith 
6029b94acceSBarry Smith    Output Parameter:
6039b94acceSBarry Smith .  sles - the SLES context
6049b94acceSBarry Smith 
6059b94acceSBarry Smith    Notes:
6069b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
6079b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
6089b94acceSBarry Smith    KSP and PC contexts as well.
6099b94acceSBarry Smith 
61036851e7fSLois Curfman McInnes    Level: beginner
61136851e7fSLois Curfman McInnes 
6129b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
6139b94acceSBarry Smith 
6149b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
6159b94acceSBarry Smith @*/
6169b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
6179b94acceSBarry Smith {
6183a40ed3dSBarry Smith   PetscFunctionBegin;
61977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
6209b94acceSBarry Smith   *sles = snes->sles;
6213a40ed3dSBarry Smith   PetscFunctionReturn(0);
6229b94acceSBarry Smith }
62382bf6240SBarry Smith 
6244a2ae208SSatish Balay #undef __FUNCT__
6254a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
626454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj)
627e24b481bSBarry Smith {
628aa482453SBarry Smith #if defined(PETSC_HAVE_AMS)
629454a90a3SBarry Smith   SNES          v = (SNES) obj;
630e24b481bSBarry Smith   int          ierr;
63143d6d2cbSBarry Smith #endif
632e24b481bSBarry Smith 
633e24b481bSBarry Smith   PetscFunctionBegin;
634e24b481bSBarry Smith 
63543d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS)
636e24b481bSBarry Smith   /* if it is already published then return */
637e24b481bSBarry Smith   if (v->amem >=0) PetscFunctionReturn(0);
638e24b481bSBarry Smith 
639454a90a3SBarry Smith   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
640e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
641e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
642e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
643e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
644454a90a3SBarry Smith   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
645433994e6SBarry Smith #endif
646e24b481bSBarry Smith   PetscFunctionReturn(0);
647e24b481bSBarry Smith }
648e24b481bSBarry Smith 
6499b94acceSBarry Smith /* -----------------------------------------------------------*/
6504a2ae208SSatish Balay #undef __FUNCT__
6514a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
6529b94acceSBarry Smith /*@C
6539b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
6549b94acceSBarry Smith 
655c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
656c7afd0dbSLois Curfman McInnes 
657c7afd0dbSLois Curfman McInnes    Input Parameters:
658c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
659c7afd0dbSLois Curfman McInnes -  type - type of method, either
660c7afd0dbSLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
661c7afd0dbSLois Curfman McInnes    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
6629b94acceSBarry Smith 
6639b94acceSBarry Smith    Output Parameter:
6649b94acceSBarry Smith .  outsnes - the new SNES context
6659b94acceSBarry Smith 
666c7afd0dbSLois Curfman McInnes    Options Database Keys:
667c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
668c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
669c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
670c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
671c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
672c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
673c1f60f51SBarry Smith 
67436851e7fSLois Curfman McInnes    Level: beginner
67536851e7fSLois Curfman McInnes 
6769b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
6779b94acceSBarry Smith 
678435da068SBarry Smith .seealso: SNESSolve(), SNESDestroy(), SNESProblemType, SNES
6799b94acceSBarry Smith @*/
6804b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
6819b94acceSBarry Smith {
6829b94acceSBarry Smith   int                 ierr;
6839b94acceSBarry Smith   SNES                snes;
6849b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
68537fcc0dbSBarry Smith 
6863a40ed3dSBarry Smith   PetscFunctionBegin;
6878ba1e511SMatthew Knepley   PetscValidPointer(outsnes);
6888ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
6898ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
6908ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
6918ba1e511SMatthew Knepley #endif
6928ba1e511SMatthew Knepley 
693d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
69429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"incorrect method type");
695d64ed03dSBarry Smith   }
6963f1db9ecSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
697b0a32e0cSBarry Smith   PetscLogObjectCreate(snes);
698e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
6999b94acceSBarry Smith   snes->max_its           = 50;
7009750a799SBarry Smith   snes->max_funcs	  = 10000;
7019b94acceSBarry Smith   snes->norm		  = 0.0;
7025a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
703b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
704b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
7059b94acceSBarry Smith     snes->atol		  = 1.e-10;
706d64ed03dSBarry Smith   } else {
707b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
708b4874afaSBarry Smith     snes->ttol            = 0.0;
709b4874afaSBarry Smith     snes->atol		  = 1.e-50;
710b4874afaSBarry Smith   }
7119b94acceSBarry Smith   snes->xtol		  = 1.e-8;
7129b94acceSBarry Smith   snes->nfuncs            = 0;
71350ffb88aSMatthew Knepley   snes->numFailures       = 0;
71450ffb88aSMatthew Knepley   snes->maxFailures       = 1;
7157a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
716639f9d9dSBarry Smith   snes->numbermonitors    = 0;
7179b94acceSBarry Smith   snes->data              = 0;
7189b94acceSBarry Smith   snes->view              = 0;
7199b94acceSBarry Smith   snes->computeumfunction = 0;
7209b94acceSBarry Smith   snes->umfunP            = 0;
7219b94acceSBarry Smith   snes->fc                = 0;
7229b94acceSBarry Smith   snes->deltatol          = 1.e-12;
7239b94acceSBarry Smith   snes->fmin              = -1.e30;
7249b94acceSBarry Smith   snes->method_class      = type;
7259b94acceSBarry Smith   snes->set_method_called = 0;
72682bf6240SBarry Smith   snes->setupcalled       = 0;
727186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
7286f24a144SLois Curfman McInnes   snes->vwork             = 0;
7296f24a144SLois Curfman McInnes   snes->nwork             = 0;
730758f92a0SBarry Smith   snes->conv_hist_len     = 0;
731758f92a0SBarry Smith   snes->conv_hist_max     = 0;
732758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
733758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
734758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
735184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
7369b94acceSBarry Smith 
7379b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
738b0a32e0cSBarry Smith   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
739b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
7409b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
7419b94acceSBarry Smith   kctx->version     = 2;
7429b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
7439b94acceSBarry Smith                              this was too large for some test cases */
7449b94acceSBarry Smith   kctx->rtol_last   = 0;
7459b94acceSBarry Smith   kctx->rtol_max    = .9;
7469b94acceSBarry Smith   kctx->gamma       = 1.0;
7479b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
7489b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
7499b94acceSBarry Smith   kctx->threshold   = .1;
7509b94acceSBarry Smith   kctx->lresid_last = 0;
7519b94acceSBarry Smith   kctx->norm_last   = 0;
7529b94acceSBarry Smith 
7539b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr);
754b0a32e0cSBarry Smith   PetscLogObjectParent(snes,snes->sles)
7555334005bSBarry Smith 
7569b94acceSBarry Smith   *outsnes = snes;
75700036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
7583a40ed3dSBarry Smith   PetscFunctionReturn(0);
7599b94acceSBarry Smith }
7609b94acceSBarry Smith 
7619b94acceSBarry Smith /* --------------------------------------------------------------- */
7624a2ae208SSatish Balay #undef __FUNCT__
7634a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
7649b94acceSBarry Smith /*@C
7659b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
7669b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
7679b94acceSBarry Smith    equations.
7689b94acceSBarry Smith 
769fee21e36SBarry Smith    Collective on SNES
770fee21e36SBarry Smith 
771c7afd0dbSLois Curfman McInnes    Input Parameters:
772c7afd0dbSLois Curfman McInnes +  snes - the SNES context
773c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
774c7afd0dbSLois Curfman McInnes .  r - vector to store function value
775c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
776c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7779b94acceSBarry Smith 
778c7afd0dbSLois Curfman McInnes    Calling sequence of func:
7798d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
780c7afd0dbSLois Curfman McInnes 
781313e4042SLois Curfman McInnes .  f - function vector
782c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
7839b94acceSBarry Smith 
7849b94acceSBarry Smith    Notes:
7859b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
7869b94acceSBarry Smith $      f'(x) x = -f(x),
787c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
7889b94acceSBarry Smith 
7899b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
7909b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
7919b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
7929b94acceSBarry Smith 
79336851e7fSLois Curfman McInnes    Level: beginner
79436851e7fSLois Curfman McInnes 
7959b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7969b94acceSBarry Smith 
797a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
7989b94acceSBarry Smith @*/
7995334005bSBarry Smith int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8009b94acceSBarry Smith {
8013a40ed3dSBarry Smith   PetscFunctionBegin;
80277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
803184914b5SBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE);
804184914b5SBarry Smith   PetscCheckSameComm(snes,r);
805a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
80629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
807a8c6a408SBarry Smith   }
808184914b5SBarry Smith 
8099b94acceSBarry Smith   snes->computefunction     = func;
8109b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8119b94acceSBarry Smith   snes->funP                = ctx;
8123a40ed3dSBarry Smith   PetscFunctionReturn(0);
8139b94acceSBarry Smith }
8149b94acceSBarry Smith 
8154a2ae208SSatish Balay #undef __FUNCT__
8164a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
8179b94acceSBarry Smith /*@
81836851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
8199b94acceSBarry Smith                          SNESSetFunction().
8209b94acceSBarry Smith 
821c7afd0dbSLois Curfman McInnes    Collective on SNES
822c7afd0dbSLois Curfman McInnes 
8239b94acceSBarry Smith    Input Parameters:
824c7afd0dbSLois Curfman McInnes +  snes - the SNES context
825c7afd0dbSLois Curfman McInnes -  x - input vector
8269b94acceSBarry Smith 
8279b94acceSBarry Smith    Output Parameter:
8283638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
8299b94acceSBarry Smith 
8301bffabb2SLois Curfman McInnes    Notes:
8319b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
8329b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
8339b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
8349b94acceSBarry Smith 
83536851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
83636851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
83736851e7fSLois Curfman McInnes    themselves.
83836851e7fSLois Curfman McInnes 
83936851e7fSLois Curfman McInnes    Level: developer
84036851e7fSLois Curfman McInnes 
8419b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
8429b94acceSBarry Smith 
843a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
8449b94acceSBarry Smith @*/
8459b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x,Vec y)
8469b94acceSBarry Smith {
8479b94acceSBarry Smith   int    ierr;
8489b94acceSBarry Smith 
8493a40ed3dSBarry Smith   PetscFunctionBegin;
850184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
851184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
852184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
853184914b5SBarry Smith   PetscCheckSameComm(snes,x);
854184914b5SBarry Smith   PetscCheckSameComm(snes,y);
855d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
85629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
857d4bb536fSBarry Smith   }
858184914b5SBarry Smith 
859d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
860d64ed03dSBarry Smith   PetscStackPush("SNES user function");
8619b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
862d64ed03dSBarry Smith   PetscStackPop;
863ae3c334cSLois Curfman McInnes   snes->nfuncs++;
864d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8653a40ed3dSBarry Smith   PetscFunctionReturn(0);
8669b94acceSBarry Smith }
8679b94acceSBarry Smith 
8684a2ae208SSatish Balay #undef __FUNCT__
8694a2ae208SSatish Balay #define __FUNCT__ "SNESSetMinimizationFunction"
8709b94acceSBarry Smith /*@C
8719b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
8729b94acceSBarry Smith    unconstrained minimization.
8739b94acceSBarry Smith 
874fee21e36SBarry Smith    Collective on SNES
875fee21e36SBarry Smith 
876c7afd0dbSLois Curfman McInnes    Input Parameters:
877c7afd0dbSLois Curfman McInnes +  snes - the SNES context
878c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
879c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
880c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
8819b94acceSBarry Smith 
882c7afd0dbSLois Curfman McInnes    Calling sequence of func:
883329f5518SBarry Smith $     func (SNES snes,Vec x,PetscReal *f,void *ctx);
884c7afd0dbSLois Curfman McInnes 
885c7afd0dbSLois Curfman McInnes +  x - input vector
8869b94acceSBarry Smith .  f - function
887c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined function context
8889b94acceSBarry Smith 
88936851e7fSLois Curfman McInnes    Level: beginner
89036851e7fSLois Curfman McInnes 
8919b94acceSBarry Smith    Notes:
8929b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
8939b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
8949b94acceSBarry Smith    SNESSetFunction().
8959b94acceSBarry Smith 
8969b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
8979b94acceSBarry Smith 
898a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
899a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
9009b94acceSBarry Smith @*/
901329f5518SBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,PetscReal*,void*),void *ctx)
9029b94acceSBarry Smith {
9033a40ed3dSBarry Smith   PetscFunctionBegin;
90477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
905a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
90629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
907a8c6a408SBarry Smith   }
9089b94acceSBarry Smith   snes->computeumfunction   = func;
9099b94acceSBarry Smith   snes->umfunP              = ctx;
9103a40ed3dSBarry Smith   PetscFunctionReturn(0);
9119b94acceSBarry Smith }
9129b94acceSBarry Smith 
9134a2ae208SSatish Balay #undef __FUNCT__
9144a2ae208SSatish Balay #define __FUNCT__ "SNESComputeMinimizationFunction"
9159b94acceSBarry Smith /*@
9169b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
9179b94acceSBarry Smith    set with SNESSetMinimizationFunction().
9189b94acceSBarry Smith 
919c7afd0dbSLois Curfman McInnes    Collective on SNES
920c7afd0dbSLois Curfman McInnes 
9219b94acceSBarry Smith    Input Parameters:
922c7afd0dbSLois Curfman McInnes +  snes - the SNES context
923c7afd0dbSLois Curfman McInnes -  x - input vector
9249b94acceSBarry Smith 
9259b94acceSBarry Smith    Output Parameter:
9269b94acceSBarry Smith .  y - function value
9279b94acceSBarry Smith 
9289b94acceSBarry Smith    Notes:
9299b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
9309b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
9319b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
932a86d99e1SLois Curfman McInnes 
93336851e7fSLois Curfman McInnes    SNESComputeMinimizationFunction() is typically used within minimization
93436851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
93536851e7fSLois Curfman McInnes    themselves.
93636851e7fSLois Curfman McInnes 
93736851e7fSLois Curfman McInnes    Level: developer
93836851e7fSLois Curfman McInnes 
939a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
940a86d99e1SLois Curfman McInnes 
941a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
942a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
9439b94acceSBarry Smith @*/
944329f5518SBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,PetscReal *y)
9459b94acceSBarry Smith {
9469b94acceSBarry Smith   int    ierr;
9473a40ed3dSBarry Smith 
9483a40ed3dSBarry Smith   PetscFunctionBegin;
949184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
950184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
951184914b5SBarry Smith   PetscCheckSameComm(snes,x);
952a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
95329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
954a8c6a408SBarry Smith   }
955184914b5SBarry Smith 
956d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
957d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
9589b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr);
959d64ed03dSBarry Smith   PetscStackPop;
960ae3c334cSLois Curfman McInnes   snes->nfuncs++;
961d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
9623a40ed3dSBarry Smith   PetscFunctionReturn(0);
9639b94acceSBarry Smith }
9649b94acceSBarry Smith 
9654a2ae208SSatish Balay #undef __FUNCT__
9664a2ae208SSatish Balay #define __FUNCT__ "SNESSetGradient"
9679b94acceSBarry Smith /*@C
9689b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
9699b94acceSBarry Smith    vector for use by the SNES routines.
9709b94acceSBarry Smith 
971c7afd0dbSLois Curfman McInnes    Collective on SNES
972c7afd0dbSLois Curfman McInnes 
9739b94acceSBarry Smith    Input Parameters:
974c7afd0dbSLois Curfman McInnes +  snes - the SNES context
9759b94acceSBarry Smith .  func - function evaluation routine
976044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
977044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
978c7afd0dbSLois Curfman McInnes -  r - vector to store gradient value
979fee21e36SBarry Smith 
9809b94acceSBarry Smith    Calling sequence of func:
9818d76a1e5SLois Curfman McInnes $     func (SNES, Vec x, Vec g, void *ctx);
9829b94acceSBarry Smith 
983c7afd0dbSLois Curfman McInnes +  x - input vector
9849b94acceSBarry Smith .  g - gradient vector
985c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined gradient context
9869b94acceSBarry Smith 
9879b94acceSBarry Smith    Notes:
9889b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
9899b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
9909b94acceSBarry Smith    SNESSetFunction().
9919b94acceSBarry Smith 
99236851e7fSLois Curfman McInnes    Level: beginner
99336851e7fSLois Curfman McInnes 
9949b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
9959b94acceSBarry Smith 
996a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
997a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
9989b94acceSBarry Smith @*/
99974679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
10009b94acceSBarry Smith {
10013a40ed3dSBarry Smith   PetscFunctionBegin;
100277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1003184914b5SBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE);
1004184914b5SBarry Smith   PetscCheckSameComm(snes,r);
1005a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
100629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1007a8c6a408SBarry Smith   }
10089b94acceSBarry Smith   snes->computefunction     = func;
10099b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
10109b94acceSBarry Smith   snes->funP                = ctx;
10113a40ed3dSBarry Smith   PetscFunctionReturn(0);
10129b94acceSBarry Smith }
10139b94acceSBarry Smith 
10144a2ae208SSatish Balay #undef __FUNCT__
10154a2ae208SSatish Balay #define __FUNCT__ "SNESComputeGradient"
10169b94acceSBarry Smith /*@
1017a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
1018a86d99e1SLois Curfman McInnes    SNESSetGradient().
10199b94acceSBarry Smith 
1020c7afd0dbSLois Curfman McInnes    Collective on SNES
1021c7afd0dbSLois Curfman McInnes 
10229b94acceSBarry Smith    Input Parameters:
1023c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1024c7afd0dbSLois Curfman McInnes -  x - input vector
10259b94acceSBarry Smith 
10269b94acceSBarry Smith    Output Parameter:
10279b94acceSBarry Smith .  y - gradient vector
10289b94acceSBarry Smith 
10299b94acceSBarry Smith    Notes:
10309b94acceSBarry Smith    SNESComputeGradient() is valid only for
10319b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
10329b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
1033a86d99e1SLois Curfman McInnes 
103436851e7fSLois Curfman McInnes    SNESComputeGradient() is typically used within minimization
103536851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
103636851e7fSLois Curfman McInnes    themselves.
103736851e7fSLois Curfman McInnes 
103836851e7fSLois Curfman McInnes    Level: developer
103936851e7fSLois Curfman McInnes 
1040a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
1041a86d99e1SLois Curfman McInnes 
1042a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
1043a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
10449b94acceSBarry Smith @*/
10459b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x,Vec y)
10469b94acceSBarry Smith {
10479b94acceSBarry Smith   int    ierr;
10483a40ed3dSBarry Smith 
10493a40ed3dSBarry Smith   PetscFunctionBegin;
1050184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1051184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1052184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
1053184914b5SBarry Smith   PetscCheckSameComm(snes,x);
1054184914b5SBarry Smith   PetscCheckSameComm(snes,y);
10553a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
105629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10573a40ed3dSBarry Smith   }
1058d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
1059d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
10609b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
1061d64ed03dSBarry Smith   PetscStackPop;
1062d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
10633a40ed3dSBarry Smith   PetscFunctionReturn(0);
10649b94acceSBarry Smith }
10659b94acceSBarry Smith 
10664a2ae208SSatish Balay #undef __FUNCT__
10674a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
106862fef451SLois Curfman McInnes /*@
106962fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
107062fef451SLois Curfman McInnes    set with SNESSetJacobian().
107162fef451SLois Curfman McInnes 
1072c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1073c7afd0dbSLois Curfman McInnes 
107462fef451SLois Curfman McInnes    Input Parameters:
1075c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1076c7afd0dbSLois Curfman McInnes -  x - input vector
107762fef451SLois Curfman McInnes 
107862fef451SLois Curfman McInnes    Output Parameters:
1079c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
108062fef451SLois Curfman McInnes .  B - optional preconditioning matrix
1081c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
1082fee21e36SBarry Smith 
108362fef451SLois Curfman McInnes    Notes:
108462fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
108562fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
108662fef451SLois Curfman McInnes 
1087dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
1088dc5a77f8SLois Curfman McInnes    flag parameter.
108962fef451SLois Curfman McInnes 
109062fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
109162fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
1092005c665bSBarry Smith    methods is SNESComputeHessian().
109362fef451SLois Curfman McInnes 
109436851e7fSLois Curfman McInnes    Level: developer
109536851e7fSLois Curfman McInnes 
109662fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
109762fef451SLois Curfman McInnes 
109862fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
109962fef451SLois Curfman McInnes @*/
11009b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
11019b94acceSBarry Smith {
11029b94acceSBarry Smith   int    ierr;
11033a40ed3dSBarry Smith 
11043a40ed3dSBarry Smith   PetscFunctionBegin;
1105184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1106184914b5SBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
1107184914b5SBarry Smith   PetscCheckSameComm(snes,X);
11083a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
110929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
11103a40ed3dSBarry Smith   }
11113a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
1112d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1113c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
1114d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
11159b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1116d64ed03dSBarry Smith   PetscStackPop;
1117d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
11186d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
111977c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
112077c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
11213a40ed3dSBarry Smith   PetscFunctionReturn(0);
11229b94acceSBarry Smith }
11239b94acceSBarry Smith 
11244a2ae208SSatish Balay #undef __FUNCT__
11254a2ae208SSatish Balay #define __FUNCT__ "SNESComputeHessian"
112662fef451SLois Curfman McInnes /*@
112762fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
112862fef451SLois Curfman McInnes    set with SNESSetHessian().
112962fef451SLois Curfman McInnes 
1130c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1131c7afd0dbSLois Curfman McInnes 
113262fef451SLois Curfman McInnes    Input Parameters:
1133c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1134c7afd0dbSLois Curfman McInnes -  x - input vector
113562fef451SLois Curfman McInnes 
113662fef451SLois Curfman McInnes    Output Parameters:
1137c7afd0dbSLois Curfman McInnes +  A - Hessian matrix
113862fef451SLois Curfman McInnes .  B - optional preconditioning matrix
1139c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
1140fee21e36SBarry Smith 
114162fef451SLois Curfman McInnes    Notes:
114262fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
114362fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
114462fef451SLois Curfman McInnes 
1145dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
1146dc5a77f8SLois Curfman McInnes    flag parameter.
114762fef451SLois Curfman McInnes 
114862fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
114962fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
115062fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
115162fef451SLois Curfman McInnes 
115236851e7fSLois Curfman McInnes    SNESComputeHessian() is typically used within minimization
115336851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
115436851e7fSLois Curfman McInnes    themselves.
115536851e7fSLois Curfman McInnes 
115636851e7fSLois Curfman McInnes    Level: developer
115736851e7fSLois Curfman McInnes 
115862fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
115962fef451SLois Curfman McInnes 
1160a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
1161a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
116262fef451SLois Curfman McInnes @*/
116362fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
11649b94acceSBarry Smith {
11659b94acceSBarry Smith   int    ierr;
11663a40ed3dSBarry Smith 
11673a40ed3dSBarry Smith   PetscFunctionBegin;
1168184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1169184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1170184914b5SBarry Smith   PetscCheckSameComm(snes,x);
11713a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
117229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
11733a40ed3dSBarry Smith   }
11743a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
1175d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
11760f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
1177d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
117862fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr);
1179d64ed03dSBarry Smith   PetscStackPop;
1180d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
11810f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
118277c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
118377c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
11843a40ed3dSBarry Smith   PetscFunctionReturn(0);
11859b94acceSBarry Smith }
11869b94acceSBarry Smith 
11874a2ae208SSatish Balay #undef __FUNCT__
11884a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
11899b94acceSBarry Smith /*@C
11909b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1191044dda88SLois Curfman McInnes    location to store the matrix.
11929b94acceSBarry Smith 
1193c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1194c7afd0dbSLois Curfman McInnes 
11959b94acceSBarry Smith    Input Parameters:
1196c7afd0dbSLois Curfman McInnes +  snes - the SNES context
11979b94acceSBarry Smith .  A - Jacobian matrix
11989b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
11999b94acceSBarry Smith .  func - Jacobian evaluation routine
1200c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
12012cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
12029b94acceSBarry Smith 
12039b94acceSBarry Smith    Calling sequence of func:
12048d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
12059b94acceSBarry Smith 
1206c7afd0dbSLois Curfman McInnes +  x - input vector
12079b94acceSBarry Smith .  A - Jacobian matrix
12089b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1209ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1210ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1211c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
12129b94acceSBarry Smith 
12139b94acceSBarry Smith    Notes:
1214dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
12152cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1216ac21db08SLois Curfman McInnes 
1217ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
12189b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
12199b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
12209b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
12219b94acceSBarry Smith    throughout the global iterations.
12229b94acceSBarry Smith 
122336851e7fSLois Curfman McInnes    Level: beginner
122436851e7fSLois Curfman McInnes 
12259b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
12269b94acceSBarry Smith 
1227ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
12289b94acceSBarry Smith @*/
1229454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
12309b94acceSBarry Smith {
12313a7fca6bSBarry Smith   int ierr;
12323a7fca6bSBarry Smith 
12333a40ed3dSBarry Smith   PetscFunctionBegin;
123477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
123500036973SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE);
123600036973SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE);
123700036973SBarry Smith   if (A) PetscCheckSameComm(snes,A);
123800036973SBarry Smith   if (B) PetscCheckSameComm(snes,B);
1239a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
124029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1241a8c6a408SBarry Smith   }
1242184914b5SBarry Smith 
12433a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
12443a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
12453a7fca6bSBarry Smith   if (A) {
12463a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
12479b94acceSBarry Smith     snes->jacobian = A;
12483a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
12493a7fca6bSBarry Smith   }
12503a7fca6bSBarry Smith   if (B) {
12513a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
12529b94acceSBarry Smith     snes->jacobian_pre = B;
12533a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
12543a7fca6bSBarry Smith   }
12553a40ed3dSBarry Smith   PetscFunctionReturn(0);
12569b94acceSBarry Smith }
125762fef451SLois Curfman McInnes 
12584a2ae208SSatish Balay #undef __FUNCT__
12594a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
1260c2aafc4cSSatish Balay /*@C
1261b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1262b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1263b4fd4287SBarry Smith 
1264c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1265c7afd0dbSLois Curfman McInnes 
1266b4fd4287SBarry Smith    Input Parameter:
1267b4fd4287SBarry Smith .  snes - the nonlinear solver context
1268b4fd4287SBarry Smith 
1269b4fd4287SBarry Smith    Output Parameters:
1270c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1271b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
127200036973SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
127300036973SBarry Smith -  func - location to put Jacobian function (or PETSC_NULL)
1274fee21e36SBarry Smith 
127536851e7fSLois Curfman McInnes    Level: advanced
127636851e7fSLois Curfman McInnes 
1277b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1278b4fd4287SBarry Smith @*/
127900036973SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
1280b4fd4287SBarry Smith {
12813a40ed3dSBarry Smith   PetscFunctionBegin;
1282184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12833a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
128429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
12853a40ed3dSBarry Smith   }
1286b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
1287b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
1288b4fd4287SBarry Smith   if (ctx)  *ctx  = snes->jacP;
128900036973SBarry Smith   if (func) *func = snes->computejacobian;
12903a40ed3dSBarry Smith   PetscFunctionReturn(0);
1291b4fd4287SBarry Smith }
1292b4fd4287SBarry Smith 
12934a2ae208SSatish Balay #undef __FUNCT__
12944a2ae208SSatish Balay #define __FUNCT__ "SNESSetHessian"
12959b94acceSBarry Smith /*@C
12969b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1297044dda88SLois Curfman McInnes    location to store the matrix.
12989b94acceSBarry Smith 
1299c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1300c7afd0dbSLois Curfman McInnes 
13019b94acceSBarry Smith    Input Parameters:
1302c7afd0dbSLois Curfman McInnes +  snes - the SNES context
13039b94acceSBarry Smith .  A - Hessian matrix
13049b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
13059b94acceSBarry Smith .  func - Jacobian evaluation routine
1306c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
1307313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
13089b94acceSBarry Smith 
13099b94acceSBarry Smith    Calling sequence of func:
13108d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
13119b94acceSBarry Smith 
1312c7afd0dbSLois Curfman McInnes +  x - input vector
13139b94acceSBarry Smith .  A - Hessian matrix
13149b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1315ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1316ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1317c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Hessian context
13189b94acceSBarry Smith 
13199b94acceSBarry Smith    Notes:
1320dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
13212cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1322ac21db08SLois Curfman McInnes 
13239b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
13249b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
13259b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
13269b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
13279b94acceSBarry Smith    throughout the global iterations.
13289b94acceSBarry Smith 
132936851e7fSLois Curfman McInnes    Level: beginner
133036851e7fSLois Curfman McInnes 
13319b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
13329b94acceSBarry Smith 
1333ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
13349b94acceSBarry Smith @*/
1335454a90a3SBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
13369b94acceSBarry Smith {
13373a7fca6bSBarry Smith   int ierr;
13383a7fca6bSBarry Smith 
13393a40ed3dSBarry Smith   PetscFunctionBegin;
134077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1341184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
1342184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
1343184914b5SBarry Smith   PetscCheckSameComm(snes,A);
1344184914b5SBarry Smith   PetscCheckSameComm(snes,B);
1345d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
134629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1347d4bb536fSBarry Smith   }
13483a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
13493a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
13503a7fca6bSBarry Smith   if (A) {
13513a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
13529b94acceSBarry Smith     snes->jacobian = A;
13533a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
13543a7fca6bSBarry Smith   }
13553a7fca6bSBarry Smith   if (B) {
13563a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
13579b94acceSBarry Smith     snes->jacobian_pre = B;
13583a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
13593a7fca6bSBarry Smith   }
13603a40ed3dSBarry Smith   PetscFunctionReturn(0);
13619b94acceSBarry Smith }
13629b94acceSBarry Smith 
13634a2ae208SSatish Balay #undef __FUNCT__
13644a2ae208SSatish Balay #define __FUNCT__ "SNESGetHessian"
136562fef451SLois Curfman McInnes /*@
136662fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
136762fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
136862fef451SLois Curfman McInnes 
1369c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object is parallel if SNES object is parallel
1370c7afd0dbSLois Curfman McInnes 
137162fef451SLois Curfman McInnes    Input Parameter:
137262fef451SLois Curfman McInnes .  snes - the nonlinear solver context
137362fef451SLois Curfman McInnes 
137462fef451SLois Curfman McInnes    Output Parameters:
1375c7afd0dbSLois Curfman McInnes +  A - location to stash Hessian matrix (or PETSC_NULL)
137662fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
1377c7afd0dbSLois Curfman McInnes -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1378fee21e36SBarry Smith 
137936851e7fSLois Curfman McInnes    Level: advanced
138036851e7fSLois Curfman McInnes 
138162fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
1382c7afd0dbSLois Curfman McInnes 
1383c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian
138462fef451SLois Curfman McInnes @*/
138562fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B,void **ctx)
138662fef451SLois Curfman McInnes {
13873a40ed3dSBarry Smith   PetscFunctionBegin;
1388184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13893a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
139029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
13913a40ed3dSBarry Smith   }
139262fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
139362fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
139462fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
13953a40ed3dSBarry Smith   PetscFunctionReturn(0);
139662fef451SLois Curfman McInnes }
139762fef451SLois Curfman McInnes 
13989b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
139945fc7adcSBarry Smith extern int SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
14009b94acceSBarry Smith 
14014a2ae208SSatish Balay #undef __FUNCT__
14024a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
14039b94acceSBarry Smith /*@
14049b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1405272ac6f2SLois Curfman McInnes    of a nonlinear solver.
14069b94acceSBarry Smith 
1407fee21e36SBarry Smith    Collective on SNES
1408fee21e36SBarry Smith 
1409c7afd0dbSLois Curfman McInnes    Input Parameters:
1410c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1411c7afd0dbSLois Curfman McInnes -  x - the solution vector
1412c7afd0dbSLois Curfman McInnes 
1413272ac6f2SLois Curfman McInnes    Notes:
1414272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1415272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1416272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1417272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1418272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1419272ac6f2SLois Curfman McInnes 
142036851e7fSLois Curfman McInnes    Level: advanced
142136851e7fSLois Curfman McInnes 
14229b94acceSBarry Smith .keywords: SNES, nonlinear, setup
14239b94acceSBarry Smith 
14249b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
14259b94acceSBarry Smith @*/
14268ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
14279b94acceSBarry Smith {
1428f1af5d2fSBarry Smith   int        ierr;
1429f1af5d2fSBarry Smith   PetscTruth flg;
14303a40ed3dSBarry Smith 
14313a40ed3dSBarry Smith   PetscFunctionBegin;
143277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
143377c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1434184914b5SBarry Smith   PetscCheckSameComm(snes,x);
14358ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
14369b94acceSBarry Smith 
1437b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1438c1f60f51SBarry Smith   /*
1439c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1440dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1441c1f60f51SBarry Smith   */
1442c1f60f51SBarry Smith   if (flg) {
1443c1f60f51SBarry Smith     Mat J;
14445a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
14455a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
14463a7fca6bSBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n");
14473a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
14483a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1449c1f60f51SBarry Smith   }
145045fc7adcSBarry Smith 
1451*32a4b47aSMatthew Knepley #if !defined(PETSC_USE_COMPLEX) || !defined(PETSC_USE_SINGLE)
145245fc7adcSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
145345fc7adcSBarry Smith   if (flg) {
145445fc7adcSBarry Smith     Mat J;
145545fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
145645fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
145745fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
145845fc7adcSBarry Smith   }
1459*32a4b47aSMatthew Knepley #endif
146045fc7adcSBarry Smith 
1461b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1462c1f60f51SBarry Smith   /*
1463dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1464c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1465c1f60f51SBarry Smith    */
146631615d04SLois Curfman McInnes   if (flg) {
1467272ac6f2SLois Curfman McInnes     Mat  J;
1468f3ef73ceSBarry Smith     SLES sles;
1469f3ef73ceSBarry Smith     PC   pc;
1470f3ef73ceSBarry Smith 
14715a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
14723a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
147393b98531SBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
147483e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
14753a7fca6bSBarry Smith       ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
1476d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
14773a7fca6bSBarry Smith       ierr = SNESSetHessian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
1478d4bb536fSBarry Smith     } else {
147929bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free option");
1480d4bb536fSBarry Smith     }
14813a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
14823a7fca6bSBarry Smith 
1483f3ef73ceSBarry Smith     /* force no preconditioner */
1484f3ef73ceSBarry Smith     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1485f3ef73ceSBarry Smith     ierr = SLESGetPC(sles,&pc);CHKERRQ(ierr);
1486f3ef73ceSBarry Smith     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1487272ac6f2SLois Curfman McInnes   }
1488f3ef73ceSBarry Smith 
14899b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
14906831982aSBarry Smith     PetscTruth iseqtr;
14916831982aSBarry Smith 
149229bbc08cSBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
149329bbc08cSBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
149429bbc08cSBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1495a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
149629bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1497a8c6a408SBarry Smith     }
1498a703fe33SLois Curfman McInnes 
1499a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
15006831982aSBarry Smith     ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr);
15016831982aSBarry Smith     if (snes->ksp_ewconv && !iseqtr) {
1502a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1503a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1504a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1505186905e3SBarry Smith       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1506a703fe33SLois Curfman McInnes     }
15079b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
150829bbc08cSBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
150929bbc08cSBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
1510a8c6a408SBarry Smith     if (!snes->computeumfunction) {
151129bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetMinimizationFunction() first");
1512a8c6a408SBarry Smith     }
151329bbc08cSBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetHessian()");
1514d4bb536fSBarry Smith   } else {
151529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
1516d4bb536fSBarry Smith   }
1517a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
151882bf6240SBarry Smith   snes->setupcalled = 1;
15193a40ed3dSBarry Smith   PetscFunctionReturn(0);
15209b94acceSBarry Smith }
15219b94acceSBarry Smith 
15224a2ae208SSatish Balay #undef __FUNCT__
15234a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
15249b94acceSBarry Smith /*@C
15259b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
15269b94acceSBarry Smith    with SNESCreate().
15279b94acceSBarry Smith 
1528c7afd0dbSLois Curfman McInnes    Collective on SNES
1529c7afd0dbSLois Curfman McInnes 
15309b94acceSBarry Smith    Input Parameter:
15319b94acceSBarry Smith .  snes - the SNES context
15329b94acceSBarry Smith 
153336851e7fSLois Curfman McInnes    Level: beginner
153436851e7fSLois Curfman McInnes 
15359b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
15369b94acceSBarry Smith 
153763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
15389b94acceSBarry Smith @*/
15399b94acceSBarry Smith int SNESDestroy(SNES snes)
15409b94acceSBarry Smith {
1541b8a78c4aSBarry Smith   int i,ierr;
15423a40ed3dSBarry Smith 
15433a40ed3dSBarry Smith   PetscFunctionBegin;
154477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15453a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1546d4bb536fSBarry Smith 
1547be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
15480f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1549be0abb6dSBarry Smith 
1550e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1551606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
15523a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
15533a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
15549b94acceSBarry Smith   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1555522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1556b8a78c4aSBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1557b8a78c4aSBarry Smith     if (snes->monitordestroy[i]) {
1558b8a78c4aSBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1559b8a78c4aSBarry Smith     }
1560b8a78c4aSBarry Smith   }
1561b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)snes);
15620452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
15633a40ed3dSBarry Smith   PetscFunctionReturn(0);
15649b94acceSBarry Smith }
15659b94acceSBarry Smith 
15669b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
15679b94acceSBarry Smith 
15684a2ae208SSatish Balay #undef __FUNCT__
15694a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
15709b94acceSBarry Smith /*@
1571d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
15729b94acceSBarry Smith 
1573c7afd0dbSLois Curfman McInnes    Collective on SNES
1574c7afd0dbSLois Curfman McInnes 
15759b94acceSBarry Smith    Input Parameters:
1576c7afd0dbSLois Curfman McInnes +  snes - the SNES context
157733174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
157833174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
157933174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
158033174efeSLois Curfman McInnes            of the change in the solution between steps
158133174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1582c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1583fee21e36SBarry Smith 
158433174efeSLois Curfman McInnes    Options Database Keys:
1585c7afd0dbSLois Curfman McInnes +    -snes_atol <atol> - Sets atol
1586c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1587c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1588c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1589c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
15909b94acceSBarry Smith 
1591d7a720efSLois Curfman McInnes    Notes:
15929b94acceSBarry Smith    The default maximum number of iterations is 50.
15939b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
15949b94acceSBarry Smith 
159536851e7fSLois Curfman McInnes    Level: intermediate
159636851e7fSLois Curfman McInnes 
159733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
15989b94acceSBarry Smith 
1599d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
16009b94acceSBarry Smith @*/
1601329f5518SBarry Smith int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
16029b94acceSBarry Smith {
16033a40ed3dSBarry Smith   PetscFunctionBegin;
160477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1605d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1606d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1607d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1608d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1609d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
16103a40ed3dSBarry Smith   PetscFunctionReturn(0);
16119b94acceSBarry Smith }
16129b94acceSBarry Smith 
16134a2ae208SSatish Balay #undef __FUNCT__
16144a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
16159b94acceSBarry Smith /*@
161633174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
161733174efeSLois Curfman McInnes 
1618c7afd0dbSLois Curfman McInnes    Not Collective
1619c7afd0dbSLois Curfman McInnes 
162033174efeSLois Curfman McInnes    Input Parameters:
1621c7afd0dbSLois Curfman McInnes +  snes - the SNES context
162233174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
162333174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
162433174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
162533174efeSLois Curfman McInnes            of the change in the solution between steps
162633174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1627c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1628fee21e36SBarry Smith 
162933174efeSLois Curfman McInnes    Notes:
163033174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
163133174efeSLois Curfman McInnes 
163236851e7fSLois Curfman McInnes    Level: intermediate
163336851e7fSLois Curfman McInnes 
163433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
163533174efeSLois Curfman McInnes 
163633174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
163733174efeSLois Curfman McInnes @*/
1638329f5518SBarry Smith int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
163933174efeSLois Curfman McInnes {
16403a40ed3dSBarry Smith   PetscFunctionBegin;
164133174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
164233174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
164333174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
164433174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
164533174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
164633174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
16473a40ed3dSBarry Smith   PetscFunctionReturn(0);
164833174efeSLois Curfman McInnes }
164933174efeSLois Curfman McInnes 
16504a2ae208SSatish Balay #undef __FUNCT__
16514a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
165233174efeSLois Curfman McInnes /*@
16539b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
16549b94acceSBarry Smith 
1655fee21e36SBarry Smith    Collective on SNES
1656fee21e36SBarry Smith 
1657c7afd0dbSLois Curfman McInnes    Input Parameters:
1658c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1659c7afd0dbSLois Curfman McInnes -  tol - tolerance
1660c7afd0dbSLois Curfman McInnes 
16619b94acceSBarry Smith    Options Database Key:
1662c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
16639b94acceSBarry Smith 
166436851e7fSLois Curfman McInnes    Level: intermediate
166536851e7fSLois Curfman McInnes 
16669b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
16679b94acceSBarry Smith 
1668d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
16699b94acceSBarry Smith @*/
1670329f5518SBarry Smith int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
16719b94acceSBarry Smith {
16723a40ed3dSBarry Smith   PetscFunctionBegin;
167377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16749b94acceSBarry Smith   snes->deltatol = tol;
16753a40ed3dSBarry Smith   PetscFunctionReturn(0);
16769b94acceSBarry Smith }
16779b94acceSBarry Smith 
16784a2ae208SSatish Balay #undef __FUNCT__
16794a2ae208SSatish Balay #define __FUNCT__ "SNESSetMinimizationFunctionTolerance"
16809b94acceSBarry Smith /*@
168177c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
16829b94acceSBarry Smith    for unconstrained minimization solvers.
16839b94acceSBarry Smith 
1684fee21e36SBarry Smith    Collective on SNES
1685fee21e36SBarry Smith 
1686c7afd0dbSLois Curfman McInnes    Input Parameters:
1687c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1688c7afd0dbSLois Curfman McInnes -  ftol - minimum function tolerance
1689c7afd0dbSLois Curfman McInnes 
16909b94acceSBarry Smith    Options Database Key:
1691c7afd0dbSLois Curfman McInnes .  -snes_fmin <ftol> - Sets ftol
16929b94acceSBarry Smith 
16939b94acceSBarry Smith    Note:
169477c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
16959b94acceSBarry Smith    methods only.
16969b94acceSBarry Smith 
169736851e7fSLois Curfman McInnes    Level: intermediate
169836851e7fSLois Curfman McInnes 
16999b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
17009b94acceSBarry Smith 
1701d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
17029b94acceSBarry Smith @*/
1703329f5518SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,PetscReal ftol)
17049b94acceSBarry Smith {
17053a40ed3dSBarry Smith   PetscFunctionBegin;
170677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17079b94acceSBarry Smith   snes->fmin = ftol;
17083a40ed3dSBarry Smith   PetscFunctionReturn(0);
17099b94acceSBarry Smith }
1710df9fa365SBarry Smith /*
1711df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1712df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1713df9fa365SBarry Smith    macros instead of functions
1714df9fa365SBarry Smith */
17154a2ae208SSatish Balay #undef __FUNCT__
17164a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
1717329f5518SBarry Smith int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1718ce1608b8SBarry Smith {
1719ce1608b8SBarry Smith   int ierr;
1720ce1608b8SBarry Smith 
1721ce1608b8SBarry Smith   PetscFunctionBegin;
1722184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1723ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1724ce1608b8SBarry Smith   PetscFunctionReturn(0);
1725ce1608b8SBarry Smith }
1726ce1608b8SBarry Smith 
17274a2ae208SSatish Balay #undef __FUNCT__
17284a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
1729b0a32e0cSBarry Smith int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1730df9fa365SBarry Smith {
1731df9fa365SBarry Smith   int ierr;
1732df9fa365SBarry Smith 
1733df9fa365SBarry Smith   PetscFunctionBegin;
1734df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1735df9fa365SBarry Smith   PetscFunctionReturn(0);
1736df9fa365SBarry Smith }
1737df9fa365SBarry Smith 
17384a2ae208SSatish Balay #undef __FUNCT__
17394a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
1740b0a32e0cSBarry Smith int SNESLGMonitorDestroy(PetscDrawLG draw)
1741df9fa365SBarry Smith {
1742df9fa365SBarry Smith   int ierr;
1743df9fa365SBarry Smith 
1744df9fa365SBarry Smith   PetscFunctionBegin;
1745df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1746df9fa365SBarry Smith   PetscFunctionReturn(0);
1747df9fa365SBarry Smith }
1748df9fa365SBarry Smith 
17499b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
17509b94acceSBarry Smith 
17514a2ae208SSatish Balay #undef __FUNCT__
17524a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
17539b94acceSBarry Smith /*@C
17545cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
17559b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
17569b94acceSBarry Smith    progress.
17579b94acceSBarry Smith 
1758fee21e36SBarry Smith    Collective on SNES
1759fee21e36SBarry Smith 
1760c7afd0dbSLois Curfman McInnes    Input Parameters:
1761c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1762c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1763b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1764b3006f0bSLois Curfman McInnes           monitor routine (use PETSC_NULL if no context is desitre)
1765b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1766b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
17679b94acceSBarry Smith 
1768c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1769329f5518SBarry Smith $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1770c7afd0dbSLois Curfman McInnes 
1771c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1772c7afd0dbSLois Curfman McInnes .    its - iteration number
1773c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
1774c7afd0dbSLois Curfman McInnes             for SNES_NONLINEAR_EQUATIONS methods
177540a0c1c6SLois Curfman McInnes .    norm - 2-norm gradient value (may be estimated)
1776c7afd0dbSLois Curfman McInnes             for SNES_UNCONSTRAINED_MINIMIZATION methods
177740a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
17789b94acceSBarry Smith 
17799665c990SLois Curfman McInnes    Options Database Keys:
1780c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1781c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1782c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1783c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1784c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1785c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1786c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1787c7afd0dbSLois Curfman McInnes                             the options database.
17889665c990SLois Curfman McInnes 
1789639f9d9dSBarry Smith    Notes:
17906bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
17916bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
17926bc08f3fSLois Curfman McInnes    order in which they were set.
1793639f9d9dSBarry Smith 
179436851e7fSLois Curfman McInnes    Level: intermediate
179536851e7fSLois Curfman McInnes 
17969b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
17979b94acceSBarry Smith 
17985cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
17999b94acceSBarry Smith @*/
1800329f5518SBarry Smith int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
18019b94acceSBarry Smith {
18023a40ed3dSBarry Smith   PetscFunctionBegin;
1803184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1804639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
180529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1806639f9d9dSBarry Smith   }
1807639f9d9dSBarry Smith 
1808639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1809b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1810639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
18113a40ed3dSBarry Smith   PetscFunctionReturn(0);
18129b94acceSBarry Smith }
18139b94acceSBarry Smith 
18144a2ae208SSatish Balay #undef __FUNCT__
18154a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
18165cd90555SBarry Smith /*@C
18175cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
18185cd90555SBarry Smith 
1819c7afd0dbSLois Curfman McInnes    Collective on SNES
1820c7afd0dbSLois Curfman McInnes 
18215cd90555SBarry Smith    Input Parameters:
18225cd90555SBarry Smith .  snes - the SNES context
18235cd90555SBarry Smith 
18245cd90555SBarry Smith    Options Database:
1825c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1826c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1827c7afd0dbSLois Curfman McInnes     set via the options database
18285cd90555SBarry Smith 
18295cd90555SBarry Smith    Notes:
18305cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
18315cd90555SBarry Smith 
183236851e7fSLois Curfman McInnes    Level: intermediate
183336851e7fSLois Curfman McInnes 
18345cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
18355cd90555SBarry Smith 
18365cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
18375cd90555SBarry Smith @*/
18385cd90555SBarry Smith int SNESClearMonitor(SNES snes)
18395cd90555SBarry Smith {
18405cd90555SBarry Smith   PetscFunctionBegin;
1841184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18425cd90555SBarry Smith   snes->numbermonitors = 0;
18435cd90555SBarry Smith   PetscFunctionReturn(0);
18445cd90555SBarry Smith }
18455cd90555SBarry Smith 
18464a2ae208SSatish Balay #undef __FUNCT__
18474a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
18489b94acceSBarry Smith /*@C
18499b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
18509b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
18519b94acceSBarry Smith 
1852fee21e36SBarry Smith    Collective on SNES
1853fee21e36SBarry Smith 
1854c7afd0dbSLois Curfman McInnes    Input Parameters:
1855c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1856c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1857c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1858c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
18599b94acceSBarry Smith 
1860c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1861329f5518SBarry Smith $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1862c7afd0dbSLois Curfman McInnes 
1863c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1864c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1865184914b5SBarry Smith .    reason - reason for convergence/divergence
1866c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
1867c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1868c7afd0dbSLois Curfman McInnes .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1869c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1870c7afd0dbSLois Curfman McInnes -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
18719b94acceSBarry Smith 
187236851e7fSLois Curfman McInnes    Level: advanced
187336851e7fSLois Curfman McInnes 
18749b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
18759b94acceSBarry Smith 
187640191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
187740191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
18789b94acceSBarry Smith @*/
1879329f5518SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
18809b94acceSBarry Smith {
18813a40ed3dSBarry Smith   PetscFunctionBegin;
1882184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18839b94acceSBarry Smith   (snes)->converged = func;
18849b94acceSBarry Smith   (snes)->cnvP      = cctx;
18853a40ed3dSBarry Smith   PetscFunctionReturn(0);
18869b94acceSBarry Smith }
18879b94acceSBarry Smith 
18884a2ae208SSatish Balay #undef __FUNCT__
18894a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
1890184914b5SBarry Smith /*@C
1891184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1892184914b5SBarry Smith 
1893184914b5SBarry Smith    Not Collective
1894184914b5SBarry Smith 
1895184914b5SBarry Smith    Input Parameter:
1896184914b5SBarry Smith .  snes - the SNES context
1897184914b5SBarry Smith 
1898184914b5SBarry Smith    Output Parameter:
1899e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1900184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1901184914b5SBarry Smith 
1902184914b5SBarry Smith    Level: intermediate
1903184914b5SBarry Smith 
1904184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1905184914b5SBarry Smith 
1906184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1907184914b5SBarry Smith 
1908184914b5SBarry Smith .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1909435da068SBarry Smith           SNESConverged_UM_LS(), SNESConverged_UM_TR(), SNESConvergedReason
1910184914b5SBarry Smith @*/
1911184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1912184914b5SBarry Smith {
1913184914b5SBarry Smith   PetscFunctionBegin;
1914184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1915184914b5SBarry Smith   *reason = snes->reason;
1916184914b5SBarry Smith   PetscFunctionReturn(0);
1917184914b5SBarry Smith }
1918184914b5SBarry Smith 
19194a2ae208SSatish Balay #undef __FUNCT__
19204a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1921c9005455SLois Curfman McInnes /*@
1922c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1923c9005455SLois Curfman McInnes 
1924fee21e36SBarry Smith    Collective on SNES
1925fee21e36SBarry Smith 
1926c7afd0dbSLois Curfman McInnes    Input Parameters:
1927c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1928c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1929cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1930758f92a0SBarry Smith .  na  - size of a and its
193164731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1932758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1933c7afd0dbSLois Curfman McInnes 
1934c9005455SLois Curfman McInnes    Notes:
1935c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1936c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1937c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1938c9005455SLois Curfman McInnes    at each step.
1939c9005455SLois Curfman McInnes 
1940c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1941c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1942c9005455SLois Curfman McInnes    during the section of code that is being timed.
1943c9005455SLois Curfman McInnes 
194436851e7fSLois Curfman McInnes    Level: intermediate
194536851e7fSLois Curfman McInnes 
1946c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1947758f92a0SBarry Smith 
194808405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1949758f92a0SBarry Smith 
1950c9005455SLois Curfman McInnes @*/
1951329f5518SBarry Smith int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset)
1952c9005455SLois Curfman McInnes {
19533a40ed3dSBarry Smith   PetscFunctionBegin;
1954c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1955c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1956c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1957758f92a0SBarry Smith   snes->conv_hist_its   = its;
1958758f92a0SBarry Smith   snes->conv_hist_max   = na;
1959758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1960758f92a0SBarry Smith   PetscFunctionReturn(0);
1961758f92a0SBarry Smith }
1962758f92a0SBarry Smith 
19634a2ae208SSatish Balay #undef __FUNCT__
19644a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
19650c4c9dddSBarry Smith /*@C
1966758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1967758f92a0SBarry Smith 
1968758f92a0SBarry Smith    Collective on SNES
1969758f92a0SBarry Smith 
1970758f92a0SBarry Smith    Input Parameter:
1971758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1972758f92a0SBarry Smith 
1973758f92a0SBarry Smith    Output Parameters:
1974758f92a0SBarry Smith .  a   - array to hold history
1975758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1976758f92a0SBarry Smith          negative if not converged) for each solve.
1977758f92a0SBarry Smith -  na  - size of a and its
1978758f92a0SBarry Smith 
1979758f92a0SBarry Smith    Notes:
1980758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1981758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1982758f92a0SBarry Smith 
1983758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1984758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1985758f92a0SBarry Smith    during the section of code that is being timed.
1986758f92a0SBarry Smith 
1987758f92a0SBarry Smith    Level: intermediate
1988758f92a0SBarry Smith 
1989758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1990758f92a0SBarry Smith 
1991758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1992758f92a0SBarry Smith 
1993758f92a0SBarry Smith @*/
1994329f5518SBarry Smith int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na)
1995758f92a0SBarry Smith {
1996758f92a0SBarry Smith   PetscFunctionBegin;
1997758f92a0SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1998758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1999758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
2000758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
20013a40ed3dSBarry Smith   PetscFunctionReturn(0);
2002c9005455SLois Curfman McInnes }
2003c9005455SLois Curfman McInnes 
2004e74ef692SMatthew Knepley #undef __FUNCT__
2005e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC"
200676b2cf59SMatthew Knepley /*@
200776b2cf59SMatthew Knepley   SNESSetRhsBC - Sets the function which applies boundary conditions
200876b2cf59SMatthew Knepley   to the Rhs of each system.
200976b2cf59SMatthew Knepley 
201076b2cf59SMatthew Knepley   Collective on SNES
201176b2cf59SMatthew Knepley 
201276b2cf59SMatthew Knepley   Input Parameters:
201376b2cf59SMatthew Knepley . snes - The nonlinear solver context
201476b2cf59SMatthew Knepley . func - The function
201576b2cf59SMatthew Knepley 
201676b2cf59SMatthew Knepley   Calling sequence of func:
201776b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx);
201876b2cf59SMatthew Knepley 
201976b2cf59SMatthew Knepley . rhs - The current rhs vector
202076b2cf59SMatthew Knepley . ctx - The user-context
202176b2cf59SMatthew Knepley 
202276b2cf59SMatthew Knepley   Level: intermediate
202376b2cf59SMatthew Knepley 
202476b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
202576b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
202676b2cf59SMatthew Knepley @*/
202776b2cf59SMatthew Knepley int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *))
202876b2cf59SMatthew Knepley {
202976b2cf59SMatthew Knepley   PetscFunctionBegin;
203076b2cf59SMatthew Knepley   PetscValidHeaderSpecific(snes, SNES_COOKIE);
203176b2cf59SMatthew Knepley   snes->applyrhsbc = func;
203276b2cf59SMatthew Knepley   PetscFunctionReturn(0);
203376b2cf59SMatthew Knepley }
203476b2cf59SMatthew Knepley 
2035e74ef692SMatthew Knepley #undef __FUNCT__
2036e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC"
203776b2cf59SMatthew Knepley /*@
203876b2cf59SMatthew Knepley   SNESDefaultRhsBC - The default boundary condition function which does nothing.
203976b2cf59SMatthew Knepley 
204076b2cf59SMatthew Knepley   Not collective
204176b2cf59SMatthew Knepley 
204276b2cf59SMatthew Knepley   Input Parameters:
204376b2cf59SMatthew Knepley . snes - The nonlinear solver context
204476b2cf59SMatthew Knepley . rhs  - The Rhs
204576b2cf59SMatthew Knepley . ctx  - The user-context
204676b2cf59SMatthew Knepley 
204776b2cf59SMatthew Knepley   Level: beginner
204876b2cf59SMatthew Knepley 
204976b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
205076b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
205176b2cf59SMatthew Knepley @*/
205276b2cf59SMatthew Knepley int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
205376b2cf59SMatthew Knepley {
205476b2cf59SMatthew Knepley   PetscFunctionBegin;
205576b2cf59SMatthew Knepley   PetscFunctionReturn(0);
205676b2cf59SMatthew Knepley }
205776b2cf59SMatthew Knepley 
2058e74ef692SMatthew Knepley #undef __FUNCT__
2059e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC"
206076b2cf59SMatthew Knepley /*@
206176b2cf59SMatthew Knepley   SNESSetSolutionBC - Sets the function which applies boundary conditions
206276b2cf59SMatthew Knepley   to the solution of each system.
206376b2cf59SMatthew Knepley 
206476b2cf59SMatthew Knepley   Collective on SNES
206576b2cf59SMatthew Knepley 
206676b2cf59SMatthew Knepley   Input Parameters:
206776b2cf59SMatthew Knepley . snes - The nonlinear solver context
206876b2cf59SMatthew Knepley . func - The function
206976b2cf59SMatthew Knepley 
207076b2cf59SMatthew Knepley   Calling sequence of func:
207176b2cf59SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
207276b2cf59SMatthew Knepley 
207376b2cf59SMatthew Knepley . sol - The current solution vector
207476b2cf59SMatthew Knepley . ctx - The user-context
207576b2cf59SMatthew Knepley 
207676b2cf59SMatthew Knepley   Level: intermediate
207776b2cf59SMatthew Knepley 
207876b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
207976b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
208076b2cf59SMatthew Knepley @*/
208176b2cf59SMatthew Knepley int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *))
208276b2cf59SMatthew Knepley {
208376b2cf59SMatthew Knepley   PetscFunctionBegin;
208476b2cf59SMatthew Knepley   PetscValidHeaderSpecific(snes, SNES_COOKIE);
208576b2cf59SMatthew Knepley   snes->applysolbc = func;
208676b2cf59SMatthew Knepley   PetscFunctionReturn(0);
208776b2cf59SMatthew Knepley }
208876b2cf59SMatthew Knepley 
2089e74ef692SMatthew Knepley #undef __FUNCT__
2090e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC"
209176b2cf59SMatthew Knepley /*@
209276b2cf59SMatthew Knepley   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
209376b2cf59SMatthew Knepley 
209476b2cf59SMatthew Knepley   Not collective
209576b2cf59SMatthew Knepley 
209676b2cf59SMatthew Knepley   Input Parameters:
209776b2cf59SMatthew Knepley . snes - The nonlinear solver context
209876b2cf59SMatthew Knepley . sol  - The solution
209976b2cf59SMatthew Knepley . ctx  - The user-context
210076b2cf59SMatthew Knepley 
210176b2cf59SMatthew Knepley   Level: beginner
210276b2cf59SMatthew Knepley 
210376b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
210476b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
210576b2cf59SMatthew Knepley @*/
210676b2cf59SMatthew Knepley int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
210776b2cf59SMatthew Knepley {
210876b2cf59SMatthew Knepley   PetscFunctionBegin;
210976b2cf59SMatthew Knepley   PetscFunctionReturn(0);
211076b2cf59SMatthew Knepley }
211176b2cf59SMatthew Knepley 
2112e74ef692SMatthew Knepley #undef __FUNCT__
2113e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
211476b2cf59SMatthew Knepley /*@
211576b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
211676b2cf59SMatthew Knepley   at the beginning of every step of the iteration.
211776b2cf59SMatthew Knepley 
211876b2cf59SMatthew Knepley   Collective on SNES
211976b2cf59SMatthew Knepley 
212076b2cf59SMatthew Knepley   Input Parameters:
212176b2cf59SMatthew Knepley . snes - The nonlinear solver context
212276b2cf59SMatthew Knepley . func - The function
212376b2cf59SMatthew Knepley 
212476b2cf59SMatthew Knepley   Calling sequence of func:
212576b2cf59SMatthew Knepley . func (TS ts, int step);
212676b2cf59SMatthew Knepley 
212776b2cf59SMatthew Knepley . step - The current step of the iteration
212876b2cf59SMatthew Knepley 
212976b2cf59SMatthew Knepley   Level: intermediate
213076b2cf59SMatthew Knepley 
213176b2cf59SMatthew Knepley .keywords: SNES, update
213276b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
213376b2cf59SMatthew Knepley @*/
213476b2cf59SMatthew Knepley int SNESSetUpdate(SNES snes, int (*func)(SNES, int))
213576b2cf59SMatthew Knepley {
213676b2cf59SMatthew Knepley   PetscFunctionBegin;
213776b2cf59SMatthew Knepley   PetscValidHeaderSpecific(snes, SNES_COOKIE);
213876b2cf59SMatthew Knepley   snes->update = func;
213976b2cf59SMatthew Knepley   PetscFunctionReturn(0);
214076b2cf59SMatthew Knepley }
214176b2cf59SMatthew Knepley 
2142e74ef692SMatthew Knepley #undef __FUNCT__
2143e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
214476b2cf59SMatthew Knepley /*@
214576b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
214676b2cf59SMatthew Knepley 
214776b2cf59SMatthew Knepley   Not collective
214876b2cf59SMatthew Knepley 
214976b2cf59SMatthew Knepley   Input Parameters:
215076b2cf59SMatthew Knepley . snes - The nonlinear solver context
215176b2cf59SMatthew Knepley . step - The current step of the iteration
215276b2cf59SMatthew Knepley 
2153205452f4SMatthew Knepley   Level: intermediate
2154205452f4SMatthew Knepley 
215576b2cf59SMatthew Knepley .keywords: SNES, update
215676b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
215776b2cf59SMatthew Knepley @*/
215876b2cf59SMatthew Knepley int SNESDefaultUpdate(SNES snes, int step)
215976b2cf59SMatthew Knepley {
216076b2cf59SMatthew Knepley   PetscFunctionBegin;
216176b2cf59SMatthew Knepley   PetscFunctionReturn(0);
216276b2cf59SMatthew Knepley }
216376b2cf59SMatthew Knepley 
21644a2ae208SSatish Balay #undef __FUNCT__
21654a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
21669b94acceSBarry Smith /*
21679b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
21689b94acceSBarry Smith    positive parameter delta.
21699b94acceSBarry Smith 
21709b94acceSBarry Smith     Input Parameters:
2171c7afd0dbSLois Curfman McInnes +   snes - the SNES context
21729b94acceSBarry Smith .   y - approximate solution of linear system
21739b94acceSBarry Smith .   fnorm - 2-norm of current function
2174c7afd0dbSLois Curfman McInnes -   delta - trust region size
21759b94acceSBarry Smith 
21769b94acceSBarry Smith     Output Parameters:
2177c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
21789b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
21799b94acceSBarry Smith     region, and exceeds zero otherwise.
2180c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
21819b94acceSBarry Smith 
21829b94acceSBarry Smith     Note:
21836831982aSBarry Smith     For non-trust region methods such as SNESEQLS, the parameter delta
21849b94acceSBarry Smith     is set to be the maximum allowable step size.
21859b94acceSBarry Smith 
21869b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
21879b94acceSBarry Smith */
2188064f8208SBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
21899b94acceSBarry Smith {
2190064f8208SBarry Smith   PetscReal   nrm;
2191ea709b57SSatish Balay   PetscScalar cnorm;
21923a40ed3dSBarry Smith   int         ierr;
21933a40ed3dSBarry Smith 
21943a40ed3dSBarry Smith   PetscFunctionBegin;
2195184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2196184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
2197184914b5SBarry Smith   PetscCheckSameComm(snes,y);
2198184914b5SBarry Smith 
2199064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
2200064f8208SBarry Smith   if (nrm > *delta) {
2201064f8208SBarry Smith      nrm = *delta/nrm;
2202064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
2203064f8208SBarry Smith      cnorm = nrm;
2204329f5518SBarry Smith      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
22059b94acceSBarry Smith      *ynorm = *delta;
22069b94acceSBarry Smith   } else {
22079b94acceSBarry Smith      *gpnorm = 0.0;
2208064f8208SBarry Smith      *ynorm = nrm;
22099b94acceSBarry Smith   }
22103a40ed3dSBarry Smith   PetscFunctionReturn(0);
22119b94acceSBarry Smith }
22129b94acceSBarry Smith 
22134a2ae208SSatish Balay #undef __FUNCT__
22144a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
22159b94acceSBarry Smith /*@
22169b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
221763a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
22189b94acceSBarry Smith 
2219c7afd0dbSLois Curfman McInnes    Collective on SNES
2220c7afd0dbSLois Curfman McInnes 
2221b2002411SLois Curfman McInnes    Input Parameters:
2222c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2223c7afd0dbSLois Curfman McInnes -  x - the solution vector
22249b94acceSBarry Smith 
22259b94acceSBarry Smith    Output Parameter:
2226b2002411SLois Curfman McInnes .  its - number of iterations until termination
22279b94acceSBarry Smith 
2228b2002411SLois Curfman McInnes    Notes:
22298ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
22308ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
22318ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
22328ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
22338ddd3da0SLois Curfman McInnes 
223436851e7fSLois Curfman McInnes    Level: beginner
223536851e7fSLois Curfman McInnes 
22369b94acceSBarry Smith .keywords: SNES, nonlinear, solve
22379b94acceSBarry Smith 
223863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
22399b94acceSBarry Smith @*/
22408ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
22419b94acceSBarry Smith {
2242f1af5d2fSBarry Smith   int        ierr;
2243f1af5d2fSBarry Smith   PetscTruth flg;
2244052efed2SBarry Smith 
22453a40ed3dSBarry Smith   PetscFunctionBegin;
224677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2247184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
2248184914b5SBarry Smith   PetscCheckSameComm(snes,x);
224974679c65SBarry Smith   PetscValidIntPointer(its);
225029bbc08cSBarry Smith   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
225174637425SBarry Smith 
225282bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
2253c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
2254758f92a0SBarry Smith   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
2255d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
225650ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
22579b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
2258d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2259b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
226093b98531SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); }
22613a40ed3dSBarry Smith   PetscFunctionReturn(0);
22629b94acceSBarry Smith }
22639b94acceSBarry Smith 
22649b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
22659b94acceSBarry Smith 
22664a2ae208SSatish Balay #undef __FUNCT__
22674a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
226882bf6240SBarry Smith /*@C
22694b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
22709b94acceSBarry Smith 
2271fee21e36SBarry Smith    Collective on SNES
2272fee21e36SBarry Smith 
2273c7afd0dbSLois Curfman McInnes    Input Parameters:
2274c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2275454a90a3SBarry Smith -  type - a known method
2276c7afd0dbSLois Curfman McInnes 
2277c7afd0dbSLois Curfman McInnes    Options Database Key:
2278454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
2279c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
2280ae12b187SLois Curfman McInnes 
22819b94acceSBarry Smith    Notes:
2282e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
22836831982aSBarry Smith +    SNESEQLS - Newton's method with line search
2284c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
22856831982aSBarry Smith .    SNESEQTR - Newton's method with trust region
2286c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
22876831982aSBarry Smith .    SNESUMTR - Newton's method with trust region
2288c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
22896831982aSBarry Smith -    SNESUMLS - Newton's method with line search
2290c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
22919b94acceSBarry Smith 
2292ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
2293ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
2294ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
2295ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
2296ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
2297ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
2298ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
2299ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
2300ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
2301b0a32e0cSBarry Smith   appropriate method.
230236851e7fSLois Curfman McInnes 
230336851e7fSLois Curfman McInnes   Level: intermediate
2304a703fe33SLois Curfman McInnes 
2305454a90a3SBarry Smith .keywords: SNES, set, type
2306435da068SBarry Smith 
2307435da068SBarry Smith .seealso: SNESType, SNESCreate()
2308435da068SBarry Smith 
23099b94acceSBarry Smith @*/
2310454a90a3SBarry Smith int SNESSetType(SNES snes,SNESType type)
23119b94acceSBarry Smith {
23120f5bd95cSBarry Smith   int        ierr,(*r)(SNES);
23136831982aSBarry Smith   PetscTruth match;
23143a40ed3dSBarry Smith 
23153a40ed3dSBarry Smith   PetscFunctionBegin;
231677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
23170f5bd95cSBarry Smith   PetscValidCharPointer(type);
231882bf6240SBarry Smith 
23196831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
23200f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
232182bf6240SBarry Smith 
232282bf6240SBarry Smith   if (snes->setupcalled) {
2323e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
232482bf6240SBarry Smith     snes->data = 0;
232582bf6240SBarry Smith   }
23267d1a2b2bSBarry Smith 
23279b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
232882bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
232982bf6240SBarry Smith 
2330b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
233182bf6240SBarry Smith 
233229bbc08cSBarry Smith   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
23331548b14aSBarry Smith 
2334606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
233582bf6240SBarry Smith   snes->data = 0;
23363a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
233782bf6240SBarry Smith 
2338454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
233982bf6240SBarry Smith   snes->set_method_called = 1;
234082bf6240SBarry Smith 
23413a40ed3dSBarry Smith   PetscFunctionReturn(0);
23429b94acceSBarry Smith }
23439b94acceSBarry Smith 
2344a847f771SSatish Balay 
23459b94acceSBarry Smith /* --------------------------------------------------------------------- */
23464a2ae208SSatish Balay #undef __FUNCT__
23474a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
23489b94acceSBarry Smith /*@C
23499b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2350f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
23519b94acceSBarry Smith 
2352fee21e36SBarry Smith    Not Collective
2353fee21e36SBarry Smith 
235436851e7fSLois Curfman McInnes    Level: advanced
235536851e7fSLois Curfman McInnes 
23569b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
23579b94acceSBarry Smith 
23589b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
23599b94acceSBarry Smith @*/
2360cf256101SBarry Smith int SNESRegisterDestroy(void)
23619b94acceSBarry Smith {
236282bf6240SBarry Smith   int ierr;
236382bf6240SBarry Smith 
23643a40ed3dSBarry Smith   PetscFunctionBegin;
236582bf6240SBarry Smith   if (SNESList) {
2366b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
236782bf6240SBarry Smith     SNESList = 0;
23689b94acceSBarry Smith   }
23694c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
23703a40ed3dSBarry Smith   PetscFunctionReturn(0);
23719b94acceSBarry Smith }
23729b94acceSBarry Smith 
23734a2ae208SSatish Balay #undef __FUNCT__
23744a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
23759b94acceSBarry Smith /*@C
23769a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
23779b94acceSBarry Smith 
2378c7afd0dbSLois Curfman McInnes    Not Collective
2379c7afd0dbSLois Curfman McInnes 
23809b94acceSBarry Smith    Input Parameter:
23814b0e389bSBarry Smith .  snes - nonlinear solver context
23829b94acceSBarry Smith 
23839b94acceSBarry Smith    Output Parameter:
23843a7fca6bSBarry Smith .  type - SNES method (a character string)
23859b94acceSBarry Smith 
238636851e7fSLois Curfman McInnes    Level: intermediate
238736851e7fSLois Curfman McInnes 
2388454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
23899b94acceSBarry Smith @*/
2390454a90a3SBarry Smith int SNESGetType(SNES snes,SNESType *type)
23919b94acceSBarry Smith {
23923a40ed3dSBarry Smith   PetscFunctionBegin;
2393184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2394454a90a3SBarry Smith   *type = snes->type_name;
23953a40ed3dSBarry Smith   PetscFunctionReturn(0);
23969b94acceSBarry Smith }
23979b94acceSBarry Smith 
23984a2ae208SSatish Balay #undef __FUNCT__
23994a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
24009b94acceSBarry Smith /*@C
24019b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
24029b94acceSBarry Smith    stored.
24039b94acceSBarry Smith 
2404c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2405c7afd0dbSLois Curfman McInnes 
24069b94acceSBarry Smith    Input Parameter:
24079b94acceSBarry Smith .  snes - the SNES context
24089b94acceSBarry Smith 
24099b94acceSBarry Smith    Output Parameter:
24109b94acceSBarry Smith .  x - the solution
24119b94acceSBarry Smith 
241236851e7fSLois Curfman McInnes    Level: advanced
241336851e7fSLois Curfman McInnes 
24149b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
24159b94acceSBarry Smith 
2416a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
24179b94acceSBarry Smith @*/
24189b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
24199b94acceSBarry Smith {
24203a40ed3dSBarry Smith   PetscFunctionBegin;
242177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
24229b94acceSBarry Smith   *x = snes->vec_sol_always;
24233a40ed3dSBarry Smith   PetscFunctionReturn(0);
24249b94acceSBarry Smith }
24259b94acceSBarry Smith 
24264a2ae208SSatish Balay #undef __FUNCT__
24274a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
24289b94acceSBarry Smith /*@C
24299b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
24309b94acceSBarry Smith    stored.
24319b94acceSBarry Smith 
2432c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2433c7afd0dbSLois Curfman McInnes 
24349b94acceSBarry Smith    Input Parameter:
24359b94acceSBarry Smith .  snes - the SNES context
24369b94acceSBarry Smith 
24379b94acceSBarry Smith    Output Parameter:
24389b94acceSBarry Smith .  x - the solution update
24399b94acceSBarry Smith 
244036851e7fSLois Curfman McInnes    Level: advanced
244136851e7fSLois Curfman McInnes 
24429b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
24439b94acceSBarry Smith 
24449b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
24459b94acceSBarry Smith @*/
24469b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
24479b94acceSBarry Smith {
24483a40ed3dSBarry Smith   PetscFunctionBegin;
244977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
24509b94acceSBarry Smith   *x = snes->vec_sol_update_always;
24513a40ed3dSBarry Smith   PetscFunctionReturn(0);
24529b94acceSBarry Smith }
24539b94acceSBarry Smith 
24544a2ae208SSatish Balay #undef __FUNCT__
24554a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
24569b94acceSBarry Smith /*@C
24573638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
24589b94acceSBarry Smith 
2459c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2460c7afd0dbSLois Curfman McInnes 
24619b94acceSBarry Smith    Input Parameter:
24629b94acceSBarry Smith .  snes - the SNES context
24639b94acceSBarry Smith 
24649b94acceSBarry Smith    Output Parameter:
24657bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
246600036973SBarry Smith .  ctx - the function context (or PETSC_NULL)
246700036973SBarry Smith -  func - the function (or PETSC_NULL)
24689b94acceSBarry Smith 
24699b94acceSBarry Smith    Notes:
24709b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
24719b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
24729b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
24739b94acceSBarry Smith 
247436851e7fSLois Curfman McInnes    Level: advanced
247536851e7fSLois Curfman McInnes 
2476a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
24779b94acceSBarry Smith 
247831615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
247931615d04SLois Curfman McInnes           SNESGetGradient()
24807bf4e008SBarry Smith 
24819b94acceSBarry Smith @*/
248200036973SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
24839b94acceSBarry Smith {
24843a40ed3dSBarry Smith   PetscFunctionBegin;
248577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2486a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
248729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
2488a8c6a408SBarry Smith   }
24897bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
24907bf4e008SBarry Smith   if (ctx)  *ctx  = snes->funP;
249100036973SBarry Smith   if (func) *func = snes->computefunction;
24923a40ed3dSBarry Smith   PetscFunctionReturn(0);
24939b94acceSBarry Smith }
24949b94acceSBarry Smith 
24954a2ae208SSatish Balay #undef __FUNCT__
24964a2ae208SSatish Balay #define __FUNCT__ "SNESGetGradient"
24979b94acceSBarry Smith /*@C
24983638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
24999b94acceSBarry Smith 
2500c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2501c7afd0dbSLois Curfman McInnes 
25029b94acceSBarry Smith    Input Parameter:
25039b94acceSBarry Smith .  snes - the SNES context
25049b94acceSBarry Smith 
25059b94acceSBarry Smith    Output Parameter:
25067bf4e008SBarry Smith +  r - the gradient (or PETSC_NULL)
25077bf4e008SBarry Smith -  ctx - the gradient context (or PETSC_NULL)
25089b94acceSBarry Smith 
25099b94acceSBarry Smith    Notes:
25109b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
25119b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
25129b94acceSBarry Smith    SNESGetFunction().
25139b94acceSBarry Smith 
251436851e7fSLois Curfman McInnes    Level: advanced
251536851e7fSLois Curfman McInnes 
25169b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
25179b94acceSBarry Smith 
25187bf4e008SBarry Smith .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(),
25197bf4e008SBarry Smith           SNESSetGradient(), SNESSetFunction()
25207bf4e008SBarry Smith 
25219b94acceSBarry Smith @*/
25227bf4e008SBarry Smith int SNESGetGradient(SNES snes,Vec *r,void **ctx)
25239b94acceSBarry Smith {
25243a40ed3dSBarry Smith   PetscFunctionBegin;
252577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
25263a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
252729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
25283a40ed3dSBarry Smith   }
25297bf4e008SBarry Smith   if (r)   *r = snes->vec_func_always;
25307bf4e008SBarry Smith   if (ctx) *ctx = snes->funP;
25313a40ed3dSBarry Smith   PetscFunctionReturn(0);
25329b94acceSBarry Smith }
25339b94acceSBarry Smith 
25344a2ae208SSatish Balay #undef __FUNCT__
25354a2ae208SSatish Balay #define __FUNCT__ "SNESGetMinimizationFunction"
25367bf4e008SBarry Smith /*@C
25379b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
25389b94acceSBarry Smith    unconstrained minimization problems.
25399b94acceSBarry Smith 
2540c7afd0dbSLois Curfman McInnes    Not Collective
2541c7afd0dbSLois Curfman McInnes 
25429b94acceSBarry Smith    Input Parameter:
25439b94acceSBarry Smith .  snes - the SNES context
25449b94acceSBarry Smith 
25459b94acceSBarry Smith    Output Parameter:
25467bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
25477bf4e008SBarry Smith -  ctx - the function context (or PETSC_NULL)
25489b94acceSBarry Smith 
25499b94acceSBarry Smith    Notes:
25509b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
25519b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
25529b94acceSBarry Smith    SNESGetFunction().
25539b94acceSBarry Smith 
255436851e7fSLois Curfman McInnes    Level: advanced
255536851e7fSLois Curfman McInnes 
25569b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
25579b94acceSBarry Smith 
25587bf4e008SBarry Smith .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction()
25597bf4e008SBarry Smith 
25609b94acceSBarry Smith @*/
2561329f5518SBarry Smith int SNESGetMinimizationFunction(SNES snes,PetscReal *r,void **ctx)
25629b94acceSBarry Smith {
25633a40ed3dSBarry Smith   PetscFunctionBegin;
256477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
256574679c65SBarry Smith   PetscValidScalarPointer(r);
25663a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
256729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
25683a40ed3dSBarry Smith   }
25697bf4e008SBarry Smith   if (r)   *r = snes->fc;
25707bf4e008SBarry Smith   if (ctx) *ctx = snes->umfunP;
25713a40ed3dSBarry Smith   PetscFunctionReturn(0);
25729b94acceSBarry Smith }
25739b94acceSBarry Smith 
25744a2ae208SSatish Balay #undef __FUNCT__
25754a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
25763c7409f5SSatish Balay /*@C
25773c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2578d850072dSLois Curfman McInnes    SNES options in the database.
25793c7409f5SSatish Balay 
2580fee21e36SBarry Smith    Collective on SNES
2581fee21e36SBarry Smith 
2582c7afd0dbSLois Curfman McInnes    Input Parameter:
2583c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2584c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2585c7afd0dbSLois Curfman McInnes 
2586d850072dSLois Curfman McInnes    Notes:
2587a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2588c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2589d850072dSLois Curfman McInnes 
259036851e7fSLois Curfman McInnes    Level: advanced
259136851e7fSLois Curfman McInnes 
25923c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2593a86d99e1SLois Curfman McInnes 
2594a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
25953c7409f5SSatish Balay @*/
25963c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
25973c7409f5SSatish Balay {
25983c7409f5SSatish Balay   int ierr;
25993c7409f5SSatish Balay 
26003a40ed3dSBarry Smith   PetscFunctionBegin;
260177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2602639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
26033c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
26043a40ed3dSBarry Smith   PetscFunctionReturn(0);
26053c7409f5SSatish Balay }
26063c7409f5SSatish Balay 
26074a2ae208SSatish Balay #undef __FUNCT__
26084a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
26093c7409f5SSatish Balay /*@C
2610f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2611d850072dSLois Curfman McInnes    SNES options in the database.
26123c7409f5SSatish Balay 
2613fee21e36SBarry Smith    Collective on SNES
2614fee21e36SBarry Smith 
2615c7afd0dbSLois Curfman McInnes    Input Parameters:
2616c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2617c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2618c7afd0dbSLois Curfman McInnes 
2619d850072dSLois Curfman McInnes    Notes:
2620a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2621c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2622d850072dSLois Curfman McInnes 
262336851e7fSLois Curfman McInnes    Level: advanced
262436851e7fSLois Curfman McInnes 
26253c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2626a86d99e1SLois Curfman McInnes 
2627a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
26283c7409f5SSatish Balay @*/
26293c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
26303c7409f5SSatish Balay {
26313c7409f5SSatish Balay   int ierr;
26323c7409f5SSatish Balay 
26333a40ed3dSBarry Smith   PetscFunctionBegin;
263477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2635639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
26363c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
26373a40ed3dSBarry Smith   PetscFunctionReturn(0);
26383c7409f5SSatish Balay }
26393c7409f5SSatish Balay 
26404a2ae208SSatish Balay #undef __FUNCT__
26414a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
26429ab63eb5SSatish Balay /*@C
26433c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
26443c7409f5SSatish Balay    SNES options in the database.
26453c7409f5SSatish Balay 
2646c7afd0dbSLois Curfman McInnes    Not Collective
2647c7afd0dbSLois Curfman McInnes 
26483c7409f5SSatish Balay    Input Parameter:
26493c7409f5SSatish Balay .  snes - the SNES context
26503c7409f5SSatish Balay 
26513c7409f5SSatish Balay    Output Parameter:
26523c7409f5SSatish Balay .  prefix - pointer to the prefix string used
26533c7409f5SSatish Balay 
26549ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
26559ab63eb5SSatish Balay    sufficient length to hold the prefix.
26569ab63eb5SSatish Balay 
265736851e7fSLois Curfman McInnes    Level: advanced
265836851e7fSLois Curfman McInnes 
26593c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2660a86d99e1SLois Curfman McInnes 
2661a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
26623c7409f5SSatish Balay @*/
26633c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
26643c7409f5SSatish Balay {
26653c7409f5SSatish Balay   int ierr;
26663c7409f5SSatish Balay 
26673a40ed3dSBarry Smith   PetscFunctionBegin;
266877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2669639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
26703a40ed3dSBarry Smith   PetscFunctionReturn(0);
26713c7409f5SSatish Balay }
26723c7409f5SSatish Balay 
2673acb85ae6SSatish Balay /*MC
2674f1af5d2fSBarry Smith    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
26759b94acceSBarry Smith 
2676b2002411SLois Curfman McInnes    Synopsis:
26773a7fca6bSBarry Smith    int SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
26789b94acceSBarry Smith 
26798d76a1e5SLois Curfman McInnes    Not collective
26808d76a1e5SLois Curfman McInnes 
2681b2002411SLois Curfman McInnes    Input Parameters:
2682c7afd0dbSLois Curfman McInnes +  name_solver - name of a new user-defined solver
2683b2002411SLois Curfman McInnes .  path - path (either absolute or relative) the library containing this solver
2684b2002411SLois Curfman McInnes .  name_create - name of routine to create method context
2685c7afd0dbSLois Curfman McInnes -  routine_create - routine to create method context
26869b94acceSBarry Smith 
2687b2002411SLois Curfman McInnes    Notes:
2688f1af5d2fSBarry Smith    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
26893a40ed3dSBarry Smith 
2690b2002411SLois Curfman McInnes    If dynamic libraries are used, then the fourth input argument (routine_create)
2691b2002411SLois Curfman McInnes    is ignored.
2692b2002411SLois Curfman McInnes 
2693b9617806SBarry Smith    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
26945ba88a07SLois Curfman McInnes    and others of the form ${any_environmental_variable} occuring in pathname will be
26955ba88a07SLois Curfman McInnes    replaced with appropriate values.
26965ba88a07SLois Curfman McInnes 
2697b2002411SLois Curfman McInnes    Sample usage:
269869225d78SLois Curfman McInnes .vb
2699f1af5d2fSBarry Smith    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2700b2002411SLois Curfman McInnes                 "MySolverCreate",MySolverCreate);
270169225d78SLois Curfman McInnes .ve
2702b2002411SLois Curfman McInnes 
2703b2002411SLois Curfman McInnes    Then, your solver can be chosen with the procedural interface via
2704b2002411SLois Curfman McInnes $     SNESSetType(snes,"my_solver")
2705b2002411SLois Curfman McInnes    or at runtime via the option
2706b2002411SLois Curfman McInnes $     -snes_type my_solver
2707b2002411SLois Curfman McInnes 
270836851e7fSLois Curfman McInnes    Level: advanced
270936851e7fSLois Curfman McInnes 
2710b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register
2711b2002411SLois Curfman McInnes 
2712b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2713b2002411SLois Curfman McInnes M*/
2714b2002411SLois Curfman McInnes 
27154a2ae208SSatish Balay #undef __FUNCT__
27164a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
2717f1af5d2fSBarry Smith int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES))
2718b2002411SLois Curfman McInnes {
2719b2002411SLois Curfman McInnes   char fullname[256];
2720b2002411SLois Curfman McInnes   int  ierr;
2721b2002411SLois Curfman McInnes 
2722b2002411SLois Curfman McInnes   PetscFunctionBegin;
2723b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2724c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2725b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2726b2002411SLois Curfman McInnes }
2727