xref: /petsc/src/snes/interface/snes.c (revision d5ba7fb7d3f82433fc93946f26018f4b1c7683c8)
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;
11*d5ba7fb7SMatthew Knepley int SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval, SNES_MinimizationFunctionEval, SNES_GradientEval;
12*d5ba7fb7SMatthew 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
184b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
18582738288SBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
18682738288SBarry Smith                                solver; hence iterations will continue until max_it
1871fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
18882738288SBarry Smith                                of convergence test
18982738288SBarry Smith .  -snes_monitor - prints residual norm at each iteration
190d132466eSBarry Smith .  -snes_vecmonitor - plots solution at each iteration
191d132466eSBarry Smith .  -snes_vecmonitor_update - plots update to solution at each iteration
19282738288SBarry Smith .  -snes_xmonitor - plots residual norm at each iteration
193e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
19436851e7fSLois Curfman McInnes -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
19582738288SBarry Smith 
19682738288SBarry Smith     Options Database for Eisenstat-Walker method:
19782738288SBarry Smith +  -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence
19882738288SBarry Smith .  -snes_ksp_eq_version ver - version of  Eisenstat-Walker method
19936851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
20036851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
20136851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
20236851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
20336851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
20436851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
20582738288SBarry Smith 
20611ca99fdSLois Curfman McInnes    Notes:
20711ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
20811ca99fdSLois Curfman McInnes    the users manual.
20983e2fdc7SBarry Smith 
21036851e7fSLois Curfman McInnes    Level: beginner
21136851e7fSLois Curfman McInnes 
2129b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
2139b94acceSBarry Smith 
21469ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
2159b94acceSBarry Smith @*/
2169b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
2179b94acceSBarry Smith {
2189b94acceSBarry Smith   SLES                sles;
219186905e3SBarry Smith   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
220f1af5d2fSBarry Smith   PetscTruth          flg;
22176b2cf59SMatthew Knepley   int                 ierr, i;
222186905e3SBarry Smith   char                *deft,type[256];
2239b94acceSBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionBegin;
22577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
226ca161407SBarry Smith 
227b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
228186905e3SBarry Smith     if (snes->type_name) {
229186905e3SBarry Smith       deft = snes->type_name;
230186905e3SBarry Smith     } else {
231186905e3SBarry Smith       if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
232186905e3SBarry Smith         deft = SNESEQLS;
233186905e3SBarry Smith       } else {
234186905e3SBarry Smith         deft = SNESUMTR;
235d64ed03dSBarry Smith       }
236d64ed03dSBarry Smith     }
2374bbc92c1SBarry Smith 
238186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
239b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
240d64ed03dSBarry Smith     if (flg) {
241186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
242186905e3SBarry Smith     } else if (!snes->type_name) {
243186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
244d64ed03dSBarry Smith     }
24593c39befSBarry Smith 
24687828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
24787828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
248186905e3SBarry Smith 
24987828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
250b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
251b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
25287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_fmin","Minimization function tolerance","SNESSetMinimizationFunctionTolerance",snes->fmin,&snes->fmin,0);CHKERRQ(ierr);
253186905e3SBarry Smith 
254b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
255186905e3SBarry Smith 
256b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
25787828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
25887828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
25987828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
26087828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
26187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
26287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
263186905e3SBarry Smith 
264b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
26593c39befSBarry Smith     if (flg) {snes->converged = 0;}
266b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
2675cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
268b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
269b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
2703a7fca6bSBarry Smith     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
2713a7fca6bSBarry Smith     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
272b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
273b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
274b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
275b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
276b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2777c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
278b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
279186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
280e24b481bSBarry Smith 
281b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2823c7409f5SSatish Balay     if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
283186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
284b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
285981c4779SBarry Smith     } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
286186905e3SBarry Smith       ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeHessian,snes->funP);CHKERRQ(ierr);
287b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2889b94acceSBarry Smith     }
289639f9d9dSBarry Smith 
29076b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
29176b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);                                                             CHKERRQ(ierr);
29276b2cf59SMatthew Knepley     }
29376b2cf59SMatthew Knepley 
294186905e3SBarry Smith     if (snes->setfromoptions) {
295186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
296639f9d9dSBarry Smith     }
297639f9d9dSBarry Smith 
298b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2994bbc92c1SBarry Smith 
3009b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
3019b94acceSBarry Smith   ierr = SLESSetFromOptions(sles);CHKERRQ(ierr);
30293993e2dSLois Curfman McInnes 
3033a40ed3dSBarry Smith   PetscFunctionReturn(0);
3049b94acceSBarry Smith }
3059b94acceSBarry Smith 
306a847f771SSatish Balay 
3074a2ae208SSatish Balay #undef __FUNCT__
3084a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
3099b94acceSBarry Smith /*@
3109b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
3119b94acceSBarry Smith    the nonlinear solvers.
3129b94acceSBarry Smith 
313fee21e36SBarry Smith    Collective on SNES
314fee21e36SBarry Smith 
315c7afd0dbSLois Curfman McInnes    Input Parameters:
316c7afd0dbSLois Curfman McInnes +  snes - the SNES context
317c7afd0dbSLois Curfman McInnes -  usrP - optional user context
318c7afd0dbSLois Curfman McInnes 
31936851e7fSLois Curfman McInnes    Level: intermediate
32036851e7fSLois Curfman McInnes 
3219b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
3229b94acceSBarry Smith 
3239b94acceSBarry Smith .seealso: SNESGetApplicationContext()
3249b94acceSBarry Smith @*/
3259b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
3269b94acceSBarry Smith {
3273a40ed3dSBarry Smith   PetscFunctionBegin;
32877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3299b94acceSBarry Smith   snes->user		= usrP;
3303a40ed3dSBarry Smith   PetscFunctionReturn(0);
3319b94acceSBarry Smith }
33274679c65SBarry Smith 
3334a2ae208SSatish Balay #undef __FUNCT__
3344a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3359b94acceSBarry Smith /*@C
3369b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3379b94acceSBarry Smith    nonlinear solvers.
3389b94acceSBarry Smith 
339c7afd0dbSLois Curfman McInnes    Not Collective
340c7afd0dbSLois Curfman McInnes 
3419b94acceSBarry Smith    Input Parameter:
3429b94acceSBarry Smith .  snes - SNES context
3439b94acceSBarry Smith 
3449b94acceSBarry Smith    Output Parameter:
3459b94acceSBarry Smith .  usrP - user context
3469b94acceSBarry Smith 
34736851e7fSLois Curfman McInnes    Level: intermediate
34836851e7fSLois Curfman McInnes 
3499b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3509b94acceSBarry Smith 
3519b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3529b94acceSBarry Smith @*/
3539b94acceSBarry Smith int SNESGetApplicationContext(SNES snes,void **usrP)
3549b94acceSBarry Smith {
3553a40ed3dSBarry Smith   PetscFunctionBegin;
35677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
3579b94acceSBarry Smith   *usrP = snes->user;
3583a40ed3dSBarry Smith   PetscFunctionReturn(0);
3599b94acceSBarry Smith }
36074679c65SBarry Smith 
3614a2ae208SSatish Balay #undef __FUNCT__
3624a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3639b94acceSBarry Smith /*@
364c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
365c8228a4eSBarry Smith    at this time.
3669b94acceSBarry Smith 
367c7afd0dbSLois Curfman McInnes    Not Collective
368c7afd0dbSLois Curfman McInnes 
3699b94acceSBarry Smith    Input Parameter:
3709b94acceSBarry Smith .  snes - SNES context
3719b94acceSBarry Smith 
3729b94acceSBarry Smith    Output Parameter:
3739b94acceSBarry Smith .  iter - iteration number
3749b94acceSBarry Smith 
375c8228a4eSBarry Smith    Notes:
376c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
377c8228a4eSBarry Smith 
378c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
37908405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
38008405cd6SLois Curfman McInnes .vb
38108405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
38208405cd6SLois Curfman McInnes       if (!(it % 2)) {
38308405cd6SLois Curfman McInnes         [compute Jacobian here]
38408405cd6SLois Curfman McInnes       }
38508405cd6SLois Curfman McInnes .ve
386c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
38708405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
388c8228a4eSBarry Smith 
38936851e7fSLois Curfman McInnes    Level: intermediate
39036851e7fSLois Curfman McInnes 
3919b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3929b94acceSBarry Smith @*/
3939b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3949b94acceSBarry Smith {
3953a40ed3dSBarry Smith   PetscFunctionBegin;
39677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
39774679c65SBarry Smith   PetscValidIntPointer(iter);
3989b94acceSBarry Smith   *iter = snes->iter;
3993a40ed3dSBarry Smith   PetscFunctionReturn(0);
4009b94acceSBarry Smith }
40174679c65SBarry Smith 
4024a2ae208SSatish Balay #undef __FUNCT__
4034a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
4049b94acceSBarry Smith /*@
4059b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
4069b94acceSBarry Smith    with SNESSSetFunction().
4079b94acceSBarry Smith 
408c7afd0dbSLois Curfman McInnes    Collective on SNES
409c7afd0dbSLois Curfman McInnes 
4109b94acceSBarry Smith    Input Parameter:
4119b94acceSBarry Smith .  snes - SNES context
4129b94acceSBarry Smith 
4139b94acceSBarry Smith    Output Parameter:
4149b94acceSBarry Smith .  fnorm - 2-norm of function
4159b94acceSBarry Smith 
4169b94acceSBarry Smith    Note:
4179b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
418a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
419a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
4209b94acceSBarry Smith 
42136851e7fSLois Curfman McInnes    Level: intermediate
42236851e7fSLois Curfman McInnes 
4239b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
424a86d99e1SLois Curfman McInnes 
42508405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
4269b94acceSBarry Smith @*/
42787828ca2SBarry Smith int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
4289b94acceSBarry Smith {
4293a40ed3dSBarry Smith   PetscFunctionBegin;
43077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
43174679c65SBarry Smith   PetscValidScalarPointer(fnorm);
43274679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
43329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_NONLINEAR_EQUATIONS only");
43474679c65SBarry Smith   }
4359b94acceSBarry Smith   *fnorm = snes->norm;
4363a40ed3dSBarry Smith   PetscFunctionReturn(0);
4379b94acceSBarry Smith }
43874679c65SBarry Smith 
4394a2ae208SSatish Balay #undef __FUNCT__
4404a2ae208SSatish Balay #define __FUNCT__ "SNESGetGradientNorm"
4419b94acceSBarry Smith /*@
4429b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
4439b94acceSBarry Smith    with SNESSSetGradient().
4449b94acceSBarry Smith 
445c7afd0dbSLois Curfman McInnes    Collective on SNES
446c7afd0dbSLois Curfman McInnes 
4479b94acceSBarry Smith    Input Parameter:
4489b94acceSBarry Smith .  snes - SNES context
4499b94acceSBarry Smith 
4509b94acceSBarry Smith    Output Parameter:
4519b94acceSBarry Smith .  fnorm - 2-norm of gradient
4529b94acceSBarry Smith 
4539b94acceSBarry Smith    Note:
4549b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
455a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
456a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
4579b94acceSBarry Smith 
45836851e7fSLois Curfman McInnes    Level: intermediate
45936851e7fSLois Curfman McInnes 
4609b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
461a86d99e1SLois Curfman McInnes 
462a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
4639b94acceSBarry Smith @*/
46487828ca2SBarry Smith int SNESGetGradientNorm(SNES snes,PetscScalar *gnorm)
4659b94acceSBarry Smith {
4663a40ed3dSBarry Smith   PetscFunctionBegin;
46777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
46874679c65SBarry Smith   PetscValidScalarPointer(gnorm);
46974679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
47029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_UNCONSTRAINED_MINIMIZATION only");
47174679c65SBarry Smith   }
4729b94acceSBarry Smith   *gnorm = snes->norm;
4733a40ed3dSBarry Smith   PetscFunctionReturn(0);
4749b94acceSBarry Smith }
47574679c65SBarry Smith 
4764a2ae208SSatish Balay #undef __FUNCT__
4774a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
4789b94acceSBarry Smith /*@
4799b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4809b94acceSBarry Smith    attempted by the nonlinear solver.
4819b94acceSBarry Smith 
482c7afd0dbSLois Curfman McInnes    Not Collective
483c7afd0dbSLois Curfman McInnes 
4849b94acceSBarry Smith    Input Parameter:
4859b94acceSBarry Smith .  snes - SNES context
4869b94acceSBarry Smith 
4879b94acceSBarry Smith    Output Parameter:
4889b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4899b94acceSBarry Smith 
490c96a6f78SLois Curfman McInnes    Notes:
491c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
492c96a6f78SLois Curfman McInnes 
49336851e7fSLois Curfman McInnes    Level: intermediate
49436851e7fSLois Curfman McInnes 
4959b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4969b94acceSBarry Smith @*/
4979b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4989b94acceSBarry Smith {
4993a40ed3dSBarry Smith   PetscFunctionBegin;
50077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
50174679c65SBarry Smith   PetscValidIntPointer(nfails);
5029b94acceSBarry Smith   *nfails = snes->nfailures;
5033a40ed3dSBarry Smith   PetscFunctionReturn(0);
5049b94acceSBarry Smith }
505a847f771SSatish Balay 
5064a2ae208SSatish Balay #undef __FUNCT__
5074a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations"
508c96a6f78SLois Curfman McInnes /*@
509c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
510c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
511c96a6f78SLois Curfman McInnes 
512c7afd0dbSLois Curfman McInnes    Not Collective
513c7afd0dbSLois Curfman McInnes 
514c96a6f78SLois Curfman McInnes    Input Parameter:
515c96a6f78SLois Curfman McInnes .  snes - SNES context
516c96a6f78SLois Curfman McInnes 
517c96a6f78SLois Curfman McInnes    Output Parameter:
518c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
519c96a6f78SLois Curfman McInnes 
520c96a6f78SLois Curfman McInnes    Notes:
521c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
522c96a6f78SLois Curfman McInnes 
52336851e7fSLois Curfman McInnes    Level: intermediate
52436851e7fSLois Curfman McInnes 
525c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
526c96a6f78SLois Curfman McInnes @*/
52786bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
528c96a6f78SLois Curfman McInnes {
5293a40ed3dSBarry Smith   PetscFunctionBegin;
530c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
531c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
532c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5333a40ed3dSBarry Smith   PetscFunctionReturn(0);
534c96a6f78SLois Curfman McInnes }
535c96a6f78SLois Curfman McInnes 
5364a2ae208SSatish Balay #undef __FUNCT__
5374a2ae208SSatish Balay #define __FUNCT__ "SNESGetSLES"
5389b94acceSBarry Smith /*@C
5399b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
5409b94acceSBarry Smith 
541c7afd0dbSLois Curfman McInnes    Not Collective, but if SNES object is parallel, then SLES object is parallel
542c7afd0dbSLois Curfman McInnes 
5439b94acceSBarry Smith    Input Parameter:
5449b94acceSBarry Smith .  snes - the SNES context
5459b94acceSBarry Smith 
5469b94acceSBarry Smith    Output Parameter:
5479b94acceSBarry Smith .  sles - the SLES context
5489b94acceSBarry Smith 
5499b94acceSBarry Smith    Notes:
5509b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
5519b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5529b94acceSBarry Smith    KSP and PC contexts as well.
5539b94acceSBarry Smith 
55436851e7fSLois Curfman McInnes    Level: beginner
55536851e7fSLois Curfman McInnes 
5569b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
5579b94acceSBarry Smith 
5589b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
5599b94acceSBarry Smith @*/
5609b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
5619b94acceSBarry Smith {
5623a40ed3dSBarry Smith   PetscFunctionBegin;
56377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
5649b94acceSBarry Smith   *sles = snes->sles;
5653a40ed3dSBarry Smith   PetscFunctionReturn(0);
5669b94acceSBarry Smith }
56782bf6240SBarry Smith 
5684a2ae208SSatish Balay #undef __FUNCT__
5694a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
570454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj)
571e24b481bSBarry Smith {
572aa482453SBarry Smith #if defined(PETSC_HAVE_AMS)
573454a90a3SBarry Smith   SNES          v = (SNES) obj;
574e24b481bSBarry Smith   int          ierr;
57543d6d2cbSBarry Smith #endif
576e24b481bSBarry Smith 
577e24b481bSBarry Smith   PetscFunctionBegin;
578e24b481bSBarry Smith 
57943d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS)
580e24b481bSBarry Smith   /* if it is already published then return */
581e24b481bSBarry Smith   if (v->amem >=0) PetscFunctionReturn(0);
582e24b481bSBarry Smith 
583454a90a3SBarry Smith   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
584e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
585e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
586e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
587e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
588454a90a3SBarry Smith   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
589433994e6SBarry Smith #endif
590e24b481bSBarry Smith   PetscFunctionReturn(0);
591e24b481bSBarry Smith }
592e24b481bSBarry Smith 
5939b94acceSBarry Smith /* -----------------------------------------------------------*/
5944a2ae208SSatish Balay #undef __FUNCT__
5954a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
5969b94acceSBarry Smith /*@C
5979b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5989b94acceSBarry Smith 
599c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
600c7afd0dbSLois Curfman McInnes 
601c7afd0dbSLois Curfman McInnes    Input Parameters:
602c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
603c7afd0dbSLois Curfman McInnes -  type - type of method, either
604c7afd0dbSLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
605c7afd0dbSLois Curfman McInnes    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
6069b94acceSBarry Smith 
6079b94acceSBarry Smith    Output Parameter:
6089b94acceSBarry Smith .  outsnes - the new SNES context
6099b94acceSBarry Smith 
610c7afd0dbSLois Curfman McInnes    Options Database Keys:
611c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
612c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
613c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
614c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
615c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
616c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
617c1f60f51SBarry Smith 
61836851e7fSLois Curfman McInnes    Level: beginner
61936851e7fSLois Curfman McInnes 
6209b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
6219b94acceSBarry Smith 
622435da068SBarry Smith .seealso: SNESSolve(), SNESDestroy(), SNESProblemType, SNES
6239b94acceSBarry Smith @*/
6244b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
6259b94acceSBarry Smith {
6269b94acceSBarry Smith   int                 ierr;
6279b94acceSBarry Smith   SNES                snes;
6289b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
62937fcc0dbSBarry Smith 
6303a40ed3dSBarry Smith   PetscFunctionBegin;
6318ba1e511SMatthew Knepley   PetscValidPointer(outsnes);
6328ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
6338ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
6348ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);                                                               CHKERRQ(ierr);
6358ba1e511SMatthew Knepley #endif
6368ba1e511SMatthew Knepley 
637d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
63829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"incorrect method type");
639d64ed03dSBarry Smith   }
6403f1db9ecSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
641b0a32e0cSBarry Smith   PetscLogObjectCreate(snes);
642e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
6439b94acceSBarry Smith   snes->max_its           = 50;
6449750a799SBarry Smith   snes->max_funcs	  = 10000;
6459b94acceSBarry Smith   snes->norm		  = 0.0;
6465a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
647b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
648b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
6499b94acceSBarry Smith     snes->atol		  = 1.e-10;
650d64ed03dSBarry Smith   } else {
651b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
652b4874afaSBarry Smith     snes->ttol            = 0.0;
653b4874afaSBarry Smith     snes->atol		  = 1.e-50;
654b4874afaSBarry Smith   }
6559b94acceSBarry Smith   snes->xtol		  = 1.e-8;
6569b94acceSBarry Smith   snes->nfuncs            = 0;
6579b94acceSBarry Smith   snes->nfailures         = 0;
6587a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
659639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6609b94acceSBarry Smith   snes->data              = 0;
6619b94acceSBarry Smith   snes->view              = 0;
6629b94acceSBarry Smith   snes->computeumfunction = 0;
6639b94acceSBarry Smith   snes->umfunP            = 0;
6649b94acceSBarry Smith   snes->fc                = 0;
6659b94acceSBarry Smith   snes->deltatol          = 1.e-12;
6669b94acceSBarry Smith   snes->fmin              = -1.e30;
6679b94acceSBarry Smith   snes->method_class      = type;
6689b94acceSBarry Smith   snes->set_method_called = 0;
66982bf6240SBarry Smith   snes->setupcalled       = 0;
670186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
6716f24a144SLois Curfman McInnes   snes->vwork             = 0;
6726f24a144SLois Curfman McInnes   snes->nwork             = 0;
673758f92a0SBarry Smith   snes->conv_hist_len     = 0;
674758f92a0SBarry Smith   snes->conv_hist_max     = 0;
675758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
676758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
677758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
678184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
6799b94acceSBarry Smith 
6809b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
681b0a32e0cSBarry Smith   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
682b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
6839b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6849b94acceSBarry Smith   kctx->version     = 2;
6859b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6869b94acceSBarry Smith                              this was too large for some test cases */
6879b94acceSBarry Smith   kctx->rtol_last   = 0;
6889b94acceSBarry Smith   kctx->rtol_max    = .9;
6899b94acceSBarry Smith   kctx->gamma       = 1.0;
6909b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6919b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6929b94acceSBarry Smith   kctx->threshold   = .1;
6939b94acceSBarry Smith   kctx->lresid_last = 0;
6949b94acceSBarry Smith   kctx->norm_last   = 0;
6959b94acceSBarry Smith 
6969b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr);
697b0a32e0cSBarry Smith   PetscLogObjectParent(snes,snes->sles)
6985334005bSBarry Smith 
6999b94acceSBarry Smith   *outsnes = snes;
70000036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
7013a40ed3dSBarry Smith   PetscFunctionReturn(0);
7029b94acceSBarry Smith }
7039b94acceSBarry Smith 
7049b94acceSBarry Smith /* --------------------------------------------------------------- */
7054a2ae208SSatish Balay #undef __FUNCT__
7064a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
7079b94acceSBarry Smith /*@C
7089b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
7099b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
7109b94acceSBarry Smith    equations.
7119b94acceSBarry Smith 
712fee21e36SBarry Smith    Collective on SNES
713fee21e36SBarry Smith 
714c7afd0dbSLois Curfman McInnes    Input Parameters:
715c7afd0dbSLois Curfman McInnes +  snes - the SNES context
716c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
717c7afd0dbSLois Curfman McInnes .  r - vector to store function value
718c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
719c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7209b94acceSBarry Smith 
721c7afd0dbSLois Curfman McInnes    Calling sequence of func:
7228d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
723c7afd0dbSLois Curfman McInnes 
724313e4042SLois Curfman McInnes .  f - function vector
725c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
7269b94acceSBarry Smith 
7279b94acceSBarry Smith    Notes:
7289b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
7299b94acceSBarry Smith $      f'(x) x = -f(x),
730c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
7319b94acceSBarry Smith 
7329b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
7339b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
7349b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
7359b94acceSBarry Smith 
73636851e7fSLois Curfman McInnes    Level: beginner
73736851e7fSLois Curfman McInnes 
7389b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
7399b94acceSBarry Smith 
740a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
7419b94acceSBarry Smith @*/
7425334005bSBarry Smith int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
7439b94acceSBarry Smith {
7443a40ed3dSBarry Smith   PetscFunctionBegin;
74577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
746184914b5SBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE);
747184914b5SBarry Smith   PetscCheckSameComm(snes,r);
748a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
74929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
750a8c6a408SBarry Smith   }
751184914b5SBarry Smith 
7529b94acceSBarry Smith   snes->computefunction     = func;
7539b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7549b94acceSBarry Smith   snes->funP                = ctx;
7553a40ed3dSBarry Smith   PetscFunctionReturn(0);
7569b94acceSBarry Smith }
7579b94acceSBarry Smith 
7584a2ae208SSatish Balay #undef __FUNCT__
7594a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
7609b94acceSBarry Smith /*@
76136851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
7629b94acceSBarry Smith                          SNESSetFunction().
7639b94acceSBarry Smith 
764c7afd0dbSLois Curfman McInnes    Collective on SNES
765c7afd0dbSLois Curfman McInnes 
7669b94acceSBarry Smith    Input Parameters:
767c7afd0dbSLois Curfman McInnes +  snes - the SNES context
768c7afd0dbSLois Curfman McInnes -  x - input vector
7699b94acceSBarry Smith 
7709b94acceSBarry Smith    Output Parameter:
7713638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7729b94acceSBarry Smith 
7731bffabb2SLois Curfman McInnes    Notes:
7749b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
7759b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
7769b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
7779b94acceSBarry Smith 
77836851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
77936851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
78036851e7fSLois Curfman McInnes    themselves.
78136851e7fSLois Curfman McInnes 
78236851e7fSLois Curfman McInnes    Level: developer
78336851e7fSLois Curfman McInnes 
7849b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7859b94acceSBarry Smith 
786a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7879b94acceSBarry Smith @*/
7889b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x,Vec y)
7899b94acceSBarry Smith {
7909b94acceSBarry Smith   int    ierr;
7919b94acceSBarry Smith 
7923a40ed3dSBarry Smith   PetscFunctionBegin;
793184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
794184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
795184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
796184914b5SBarry Smith   PetscCheckSameComm(snes,x);
797184914b5SBarry Smith   PetscCheckSameComm(snes,y);
798d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
79929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
800d4bb536fSBarry Smith   }
801184914b5SBarry Smith 
802*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
803d64ed03dSBarry Smith   PetscStackPush("SNES user function");
8049b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
805d64ed03dSBarry Smith   PetscStackPop;
806ae3c334cSLois Curfman McInnes   snes->nfuncs++;
807*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
8083a40ed3dSBarry Smith   PetscFunctionReturn(0);
8099b94acceSBarry Smith }
8109b94acceSBarry Smith 
8114a2ae208SSatish Balay #undef __FUNCT__
8124a2ae208SSatish Balay #define __FUNCT__ "SNESSetMinimizationFunction"
8139b94acceSBarry Smith /*@C
8149b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
8159b94acceSBarry Smith    unconstrained minimization.
8169b94acceSBarry Smith 
817fee21e36SBarry Smith    Collective on SNES
818fee21e36SBarry Smith 
819c7afd0dbSLois Curfman McInnes    Input Parameters:
820c7afd0dbSLois Curfman McInnes +  snes - the SNES context
821c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
822c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
823c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
8249b94acceSBarry Smith 
825c7afd0dbSLois Curfman McInnes    Calling sequence of func:
826329f5518SBarry Smith $     func (SNES snes,Vec x,PetscReal *f,void *ctx);
827c7afd0dbSLois Curfman McInnes 
828c7afd0dbSLois Curfman McInnes +  x - input vector
8299b94acceSBarry Smith .  f - function
830c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined function context
8319b94acceSBarry Smith 
83236851e7fSLois Curfman McInnes    Level: beginner
83336851e7fSLois Curfman McInnes 
8349b94acceSBarry Smith    Notes:
8359b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
8369b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
8379b94acceSBarry Smith    SNESSetFunction().
8389b94acceSBarry Smith 
8399b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
8409b94acceSBarry Smith 
841a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
842a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
8439b94acceSBarry Smith @*/
844329f5518SBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,PetscReal*,void*),void *ctx)
8459b94acceSBarry Smith {
8463a40ed3dSBarry Smith   PetscFunctionBegin;
84777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
848a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
84929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
850a8c6a408SBarry Smith   }
8519b94acceSBarry Smith   snes->computeumfunction   = func;
8529b94acceSBarry Smith   snes->umfunP              = ctx;
8533a40ed3dSBarry Smith   PetscFunctionReturn(0);
8549b94acceSBarry Smith }
8559b94acceSBarry Smith 
8564a2ae208SSatish Balay #undef __FUNCT__
8574a2ae208SSatish Balay #define __FUNCT__ "SNESComputeMinimizationFunction"
8589b94acceSBarry Smith /*@
8599b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
8609b94acceSBarry Smith    set with SNESSetMinimizationFunction().
8619b94acceSBarry Smith 
862c7afd0dbSLois Curfman McInnes    Collective on SNES
863c7afd0dbSLois Curfman McInnes 
8649b94acceSBarry Smith    Input Parameters:
865c7afd0dbSLois Curfman McInnes +  snes - the SNES context
866c7afd0dbSLois Curfman McInnes -  x - input vector
8679b94acceSBarry Smith 
8689b94acceSBarry Smith    Output Parameter:
8699b94acceSBarry Smith .  y - function value
8709b94acceSBarry Smith 
8719b94acceSBarry Smith    Notes:
8729b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
8739b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8749b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
875a86d99e1SLois Curfman McInnes 
87636851e7fSLois Curfman McInnes    SNESComputeMinimizationFunction() is typically used within minimization
87736851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
87836851e7fSLois Curfman McInnes    themselves.
87936851e7fSLois Curfman McInnes 
88036851e7fSLois Curfman McInnes    Level: developer
88136851e7fSLois Curfman McInnes 
882a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
883a86d99e1SLois Curfman McInnes 
884a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
885a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
8869b94acceSBarry Smith @*/
887329f5518SBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,PetscReal *y)
8889b94acceSBarry Smith {
8899b94acceSBarry Smith   int    ierr;
8903a40ed3dSBarry Smith 
8913a40ed3dSBarry Smith   PetscFunctionBegin;
892184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
893184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
894184914b5SBarry Smith   PetscCheckSameComm(snes,x);
895a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
89629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
897a8c6a408SBarry Smith   }
898184914b5SBarry Smith 
899*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
900d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
9019b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr);
902d64ed03dSBarry Smith   PetscStackPop;
903ae3c334cSLois Curfman McInnes   snes->nfuncs++;
904*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
9053a40ed3dSBarry Smith   PetscFunctionReturn(0);
9069b94acceSBarry Smith }
9079b94acceSBarry Smith 
9084a2ae208SSatish Balay #undef __FUNCT__
9094a2ae208SSatish Balay #define __FUNCT__ "SNESSetGradient"
9109b94acceSBarry Smith /*@C
9119b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
9129b94acceSBarry Smith    vector for use by the SNES routines.
9139b94acceSBarry Smith 
914c7afd0dbSLois Curfman McInnes    Collective on SNES
915c7afd0dbSLois Curfman McInnes 
9169b94acceSBarry Smith    Input Parameters:
917c7afd0dbSLois Curfman McInnes +  snes - the SNES context
9189b94acceSBarry Smith .  func - function evaluation routine
919044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
920044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
921c7afd0dbSLois Curfman McInnes -  r - vector to store gradient value
922fee21e36SBarry Smith 
9239b94acceSBarry Smith    Calling sequence of func:
9248d76a1e5SLois Curfman McInnes $     func (SNES, Vec x, Vec g, void *ctx);
9259b94acceSBarry Smith 
926c7afd0dbSLois Curfman McInnes +  x - input vector
9279b94acceSBarry Smith .  g - gradient vector
928c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined gradient context
9299b94acceSBarry Smith 
9309b94acceSBarry Smith    Notes:
9319b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
9329b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
9339b94acceSBarry Smith    SNESSetFunction().
9349b94acceSBarry Smith 
93536851e7fSLois Curfman McInnes    Level: beginner
93636851e7fSLois Curfman McInnes 
9379b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
9389b94acceSBarry Smith 
939a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
940a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
9419b94acceSBarry Smith @*/
94274679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
9439b94acceSBarry Smith {
9443a40ed3dSBarry Smith   PetscFunctionBegin;
94577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
946184914b5SBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE);
947184914b5SBarry Smith   PetscCheckSameComm(snes,r);
948a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
94929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
950a8c6a408SBarry Smith   }
9519b94acceSBarry Smith   snes->computefunction     = func;
9529b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
9539b94acceSBarry Smith   snes->funP                = ctx;
9543a40ed3dSBarry Smith   PetscFunctionReturn(0);
9559b94acceSBarry Smith }
9569b94acceSBarry Smith 
9574a2ae208SSatish Balay #undef __FUNCT__
9584a2ae208SSatish Balay #define __FUNCT__ "SNESComputeGradient"
9599b94acceSBarry Smith /*@
960a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
961a86d99e1SLois Curfman McInnes    SNESSetGradient().
9629b94acceSBarry Smith 
963c7afd0dbSLois Curfman McInnes    Collective on SNES
964c7afd0dbSLois Curfman McInnes 
9659b94acceSBarry Smith    Input Parameters:
966c7afd0dbSLois Curfman McInnes +  snes - the SNES context
967c7afd0dbSLois Curfman McInnes -  x - input vector
9689b94acceSBarry Smith 
9699b94acceSBarry Smith    Output Parameter:
9709b94acceSBarry Smith .  y - gradient vector
9719b94acceSBarry Smith 
9729b94acceSBarry Smith    Notes:
9739b94acceSBarry Smith    SNESComputeGradient() is valid only for
9749b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
9759b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
976a86d99e1SLois Curfman McInnes 
97736851e7fSLois Curfman McInnes    SNESComputeGradient() is typically used within minimization
97836851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
97936851e7fSLois Curfman McInnes    themselves.
98036851e7fSLois Curfman McInnes 
98136851e7fSLois Curfman McInnes    Level: developer
98236851e7fSLois Curfman McInnes 
983a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
984a86d99e1SLois Curfman McInnes 
985a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
986a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
9879b94acceSBarry Smith @*/
9889b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x,Vec y)
9899b94acceSBarry Smith {
9909b94acceSBarry Smith   int    ierr;
9913a40ed3dSBarry Smith 
9923a40ed3dSBarry Smith   PetscFunctionBegin;
993184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
994184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
995184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
996184914b5SBarry Smith   PetscCheckSameComm(snes,x);
997184914b5SBarry Smith   PetscCheckSameComm(snes,y);
9983a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
99929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10003a40ed3dSBarry Smith   }
1001*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
1002d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
10039b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
1004d64ed03dSBarry Smith   PetscStackPop;
1005*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
10063a40ed3dSBarry Smith   PetscFunctionReturn(0);
10079b94acceSBarry Smith }
10089b94acceSBarry Smith 
10094a2ae208SSatish Balay #undef __FUNCT__
10104a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
101162fef451SLois Curfman McInnes /*@
101262fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
101362fef451SLois Curfman McInnes    set with SNESSetJacobian().
101462fef451SLois Curfman McInnes 
1015c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1016c7afd0dbSLois Curfman McInnes 
101762fef451SLois Curfman McInnes    Input Parameters:
1018c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1019c7afd0dbSLois Curfman McInnes -  x - input vector
102062fef451SLois Curfman McInnes 
102162fef451SLois Curfman McInnes    Output Parameters:
1022c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
102362fef451SLois Curfman McInnes .  B - optional preconditioning matrix
1024c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
1025fee21e36SBarry Smith 
102662fef451SLois Curfman McInnes    Notes:
102762fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
102862fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
102962fef451SLois Curfman McInnes 
1030dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
1031dc5a77f8SLois Curfman McInnes    flag parameter.
103262fef451SLois Curfman McInnes 
103362fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
103462fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
1035005c665bSBarry Smith    methods is SNESComputeHessian().
103662fef451SLois Curfman McInnes 
103736851e7fSLois Curfman McInnes    Level: developer
103836851e7fSLois Curfman McInnes 
103962fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
104062fef451SLois Curfman McInnes 
104162fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
104262fef451SLois Curfman McInnes @*/
10439b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
10449b94acceSBarry Smith {
10459b94acceSBarry Smith   int    ierr;
10463a40ed3dSBarry Smith 
10473a40ed3dSBarry Smith   PetscFunctionBegin;
1048184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1049184914b5SBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
1050184914b5SBarry Smith   PetscCheckSameComm(snes,X);
10513a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
105229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
10533a40ed3dSBarry Smith   }
10543a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
1055*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
1056c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
1057d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
10589b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
1059d64ed03dSBarry Smith   PetscStackPop;
1060*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
10616d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
106277c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
106377c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
10643a40ed3dSBarry Smith   PetscFunctionReturn(0);
10659b94acceSBarry Smith }
10669b94acceSBarry Smith 
10674a2ae208SSatish Balay #undef __FUNCT__
10684a2ae208SSatish Balay #define __FUNCT__ "SNESComputeHessian"
106962fef451SLois Curfman McInnes /*@
107062fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
107162fef451SLois Curfman McInnes    set with SNESSetHessian().
107262fef451SLois Curfman McInnes 
1073c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1074c7afd0dbSLois Curfman McInnes 
107562fef451SLois Curfman McInnes    Input Parameters:
1076c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1077c7afd0dbSLois Curfman McInnes -  x - input vector
107862fef451SLois Curfman McInnes 
107962fef451SLois Curfman McInnes    Output Parameters:
1080c7afd0dbSLois Curfman McInnes +  A - Hessian matrix
108162fef451SLois Curfman McInnes .  B - optional preconditioning matrix
1082c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
1083fee21e36SBarry Smith 
108462fef451SLois Curfman McInnes    Notes:
108562fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
108662fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
108762fef451SLois Curfman McInnes 
1088dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
1089dc5a77f8SLois Curfman McInnes    flag parameter.
109062fef451SLois Curfman McInnes 
109162fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
109262fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
109362fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
109462fef451SLois Curfman McInnes 
109536851e7fSLois Curfman McInnes    SNESComputeHessian() is typically used within minimization
109636851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
109736851e7fSLois Curfman McInnes    themselves.
109836851e7fSLois Curfman McInnes 
109936851e7fSLois Curfman McInnes    Level: developer
110036851e7fSLois Curfman McInnes 
110162fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
110262fef451SLois Curfman McInnes 
1103a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
1104a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
110562fef451SLois Curfman McInnes @*/
110662fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
11079b94acceSBarry Smith {
11089b94acceSBarry Smith   int    ierr;
11093a40ed3dSBarry Smith 
11103a40ed3dSBarry Smith   PetscFunctionBegin;
1111184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1112184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1113184914b5SBarry Smith   PetscCheckSameComm(snes,x);
11143a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
111529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
11163a40ed3dSBarry Smith   }
11173a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
1118*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
11190f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
1120d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
112162fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr);
1122d64ed03dSBarry Smith   PetscStackPop;
1123*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
11240f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
112577c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
112677c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
11273a40ed3dSBarry Smith   PetscFunctionReturn(0);
11289b94acceSBarry Smith }
11299b94acceSBarry Smith 
11304a2ae208SSatish Balay #undef __FUNCT__
11314a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
11329b94acceSBarry Smith /*@C
11339b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1134044dda88SLois Curfman McInnes    location to store the matrix.
11359b94acceSBarry Smith 
1136c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1137c7afd0dbSLois Curfman McInnes 
11389b94acceSBarry Smith    Input Parameters:
1139c7afd0dbSLois Curfman McInnes +  snes - the SNES context
11409b94acceSBarry Smith .  A - Jacobian matrix
11419b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
11429b94acceSBarry Smith .  func - Jacobian evaluation routine
1143c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
11442cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
11459b94acceSBarry Smith 
11469b94acceSBarry Smith    Calling sequence of func:
11478d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
11489b94acceSBarry Smith 
1149c7afd0dbSLois Curfman McInnes +  x - input vector
11509b94acceSBarry Smith .  A - Jacobian matrix
11519b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1152ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1153ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1154c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
11559b94acceSBarry Smith 
11569b94acceSBarry Smith    Notes:
1157dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
11582cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1159ac21db08SLois Curfman McInnes 
1160ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
11619b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
11629b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
11639b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
11649b94acceSBarry Smith    throughout the global iterations.
11659b94acceSBarry Smith 
116636851e7fSLois Curfman McInnes    Level: beginner
116736851e7fSLois Curfman McInnes 
11689b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
11699b94acceSBarry Smith 
1170ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
11719b94acceSBarry Smith @*/
1172454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
11739b94acceSBarry Smith {
11743a7fca6bSBarry Smith   int ierr;
11753a7fca6bSBarry Smith 
11763a40ed3dSBarry Smith   PetscFunctionBegin;
117777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
117800036973SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE);
117900036973SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE);
118000036973SBarry Smith   if (A) PetscCheckSameComm(snes,A);
118100036973SBarry Smith   if (B) PetscCheckSameComm(snes,B);
1182a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
118329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1184a8c6a408SBarry Smith   }
1185184914b5SBarry Smith 
11863a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
11873a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
11883a7fca6bSBarry Smith   if (A) {
11893a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
11909b94acceSBarry Smith     snes->jacobian = A;
11913a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
11923a7fca6bSBarry Smith   }
11933a7fca6bSBarry Smith   if (B) {
11943a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
11959b94acceSBarry Smith     snes->jacobian_pre = B;
11963a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
11973a7fca6bSBarry Smith   }
11983a40ed3dSBarry Smith   PetscFunctionReturn(0);
11999b94acceSBarry Smith }
120062fef451SLois Curfman McInnes 
12014a2ae208SSatish Balay #undef __FUNCT__
12024a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
1203c2aafc4cSSatish Balay /*@C
1204b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1205b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1206b4fd4287SBarry Smith 
1207c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1208c7afd0dbSLois Curfman McInnes 
1209b4fd4287SBarry Smith    Input Parameter:
1210b4fd4287SBarry Smith .  snes - the nonlinear solver context
1211b4fd4287SBarry Smith 
1212b4fd4287SBarry Smith    Output Parameters:
1213c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1214b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
121500036973SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
121600036973SBarry Smith -  func - location to put Jacobian function (or PETSC_NULL)
1217fee21e36SBarry Smith 
121836851e7fSLois Curfman McInnes    Level: advanced
121936851e7fSLois Curfman McInnes 
1220b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1221b4fd4287SBarry Smith @*/
122200036973SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
1223b4fd4287SBarry Smith {
12243a40ed3dSBarry Smith   PetscFunctionBegin;
1225184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12263a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
122729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
12283a40ed3dSBarry Smith   }
1229b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
1230b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
1231b4fd4287SBarry Smith   if (ctx)  *ctx  = snes->jacP;
123200036973SBarry Smith   if (func) *func = snes->computejacobian;
12333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1234b4fd4287SBarry Smith }
1235b4fd4287SBarry Smith 
12364a2ae208SSatish Balay #undef __FUNCT__
12374a2ae208SSatish Balay #define __FUNCT__ "SNESSetHessian"
12389b94acceSBarry Smith /*@C
12399b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1240044dda88SLois Curfman McInnes    location to store the matrix.
12419b94acceSBarry Smith 
1242c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1243c7afd0dbSLois Curfman McInnes 
12449b94acceSBarry Smith    Input Parameters:
1245c7afd0dbSLois Curfman McInnes +  snes - the SNES context
12469b94acceSBarry Smith .  A - Hessian matrix
12479b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
12489b94acceSBarry Smith .  func - Jacobian evaluation routine
1249c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
1250313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
12519b94acceSBarry Smith 
12529b94acceSBarry Smith    Calling sequence of func:
12538d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
12549b94acceSBarry Smith 
1255c7afd0dbSLois Curfman McInnes +  x - input vector
12569b94acceSBarry Smith .  A - Hessian matrix
12579b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1258ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1259ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1260c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Hessian context
12619b94acceSBarry Smith 
12629b94acceSBarry Smith    Notes:
1263dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
12642cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1265ac21db08SLois Curfman McInnes 
12669b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
12679b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
12689b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
12699b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
12709b94acceSBarry Smith    throughout the global iterations.
12719b94acceSBarry Smith 
127236851e7fSLois Curfman McInnes    Level: beginner
127336851e7fSLois Curfman McInnes 
12749b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
12759b94acceSBarry Smith 
1276ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
12779b94acceSBarry Smith @*/
1278454a90a3SBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
12799b94acceSBarry Smith {
12803a7fca6bSBarry Smith   int ierr;
12813a7fca6bSBarry Smith 
12823a40ed3dSBarry Smith   PetscFunctionBegin;
128377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1284184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
1285184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
1286184914b5SBarry Smith   PetscCheckSameComm(snes,A);
1287184914b5SBarry Smith   PetscCheckSameComm(snes,B);
1288d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
128929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1290d4bb536fSBarry Smith   }
12913a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
12923a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
12933a7fca6bSBarry Smith   if (A) {
12943a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
12959b94acceSBarry Smith     snes->jacobian = A;
12963a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
12973a7fca6bSBarry Smith   }
12983a7fca6bSBarry Smith   if (B) {
12993a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
13009b94acceSBarry Smith     snes->jacobian_pre = B;
13013a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
13023a7fca6bSBarry Smith   }
13033a40ed3dSBarry Smith   PetscFunctionReturn(0);
13049b94acceSBarry Smith }
13059b94acceSBarry Smith 
13064a2ae208SSatish Balay #undef __FUNCT__
13074a2ae208SSatish Balay #define __FUNCT__ "SNESGetHessian"
130862fef451SLois Curfman McInnes /*@
130962fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
131062fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
131162fef451SLois Curfman McInnes 
1312c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object is parallel if SNES object is parallel
1313c7afd0dbSLois Curfman McInnes 
131462fef451SLois Curfman McInnes    Input Parameter:
131562fef451SLois Curfman McInnes .  snes - the nonlinear solver context
131662fef451SLois Curfman McInnes 
131762fef451SLois Curfman McInnes    Output Parameters:
1318c7afd0dbSLois Curfman McInnes +  A - location to stash Hessian matrix (or PETSC_NULL)
131962fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
1320c7afd0dbSLois Curfman McInnes -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1321fee21e36SBarry Smith 
132236851e7fSLois Curfman McInnes    Level: advanced
132336851e7fSLois Curfman McInnes 
132462fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
1325c7afd0dbSLois Curfman McInnes 
1326c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian
132762fef451SLois Curfman McInnes @*/
132862fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B,void **ctx)
132962fef451SLois Curfman McInnes {
13303a40ed3dSBarry Smith   PetscFunctionBegin;
1331184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13323a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
133329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
13343a40ed3dSBarry Smith   }
133562fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
133662fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
133762fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
13383a40ed3dSBarry Smith   PetscFunctionReturn(0);
133962fef451SLois Curfman McInnes }
134062fef451SLois Curfman McInnes 
13419b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
13429b94acceSBarry Smith 
13434a2ae208SSatish Balay #undef __FUNCT__
13444a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
13459b94acceSBarry Smith /*@
13469b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1347272ac6f2SLois Curfman McInnes    of a nonlinear solver.
13489b94acceSBarry Smith 
1349fee21e36SBarry Smith    Collective on SNES
1350fee21e36SBarry Smith 
1351c7afd0dbSLois Curfman McInnes    Input Parameters:
1352c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1353c7afd0dbSLois Curfman McInnes -  x - the solution vector
1354c7afd0dbSLois Curfman McInnes 
1355272ac6f2SLois Curfman McInnes    Notes:
1356272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1357272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1358272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1359272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1360272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1361272ac6f2SLois Curfman McInnes 
136236851e7fSLois Curfman McInnes    Level: advanced
136336851e7fSLois Curfman McInnes 
13649b94acceSBarry Smith .keywords: SNES, nonlinear, setup
13659b94acceSBarry Smith 
13669b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
13679b94acceSBarry Smith @*/
13688ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
13699b94acceSBarry Smith {
1370f1af5d2fSBarry Smith   int        ierr;
1371f1af5d2fSBarry Smith   PetscTruth flg;
13723a40ed3dSBarry Smith 
13733a40ed3dSBarry Smith   PetscFunctionBegin;
137477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
137577c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1376184914b5SBarry Smith   PetscCheckSameComm(snes,x);
13778ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
13789b94acceSBarry Smith 
1379b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1380c1f60f51SBarry Smith   /*
1381c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1382dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1383c1f60f51SBarry Smith   */
1384c1f60f51SBarry Smith   if (flg) {
1385c1f60f51SBarry Smith     Mat J;
13865a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
13875a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
13883a7fca6bSBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n");
13893a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
13903a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
1391c1f60f51SBarry Smith   }
1392b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1393c1f60f51SBarry Smith   /*
1394dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1395c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1396c1f60f51SBarry Smith    */
139731615d04SLois Curfman McInnes   if (flg) {
1398272ac6f2SLois Curfman McInnes     Mat  J;
1399f3ef73ceSBarry Smith     SLES sles;
1400f3ef73ceSBarry Smith     PC   pc;
1401f3ef73ceSBarry Smith 
14025a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
14033a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
140493b98531SBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
140583e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
14063a7fca6bSBarry Smith       ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
1407d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
14083a7fca6bSBarry Smith       ierr = SNESSetHessian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
1409d4bb536fSBarry Smith     } else {
141029bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free option");
1411d4bb536fSBarry Smith     }
14123a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
14133a7fca6bSBarry Smith 
1414f3ef73ceSBarry Smith     /* force no preconditioner */
1415f3ef73ceSBarry Smith     ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1416f3ef73ceSBarry Smith     ierr = SLESGetPC(sles,&pc);CHKERRQ(ierr);
1417f3ef73ceSBarry Smith     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1418272ac6f2SLois Curfman McInnes   }
1419f3ef73ceSBarry Smith 
14209b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
14216831982aSBarry Smith     PetscTruth iseqtr;
14226831982aSBarry Smith 
142329bbc08cSBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
142429bbc08cSBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
142529bbc08cSBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1426a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
142729bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1428a8c6a408SBarry Smith     }
1429a703fe33SLois Curfman McInnes 
1430a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
14316831982aSBarry Smith     ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr);
14326831982aSBarry Smith     if (snes->ksp_ewconv && !iseqtr) {
1433a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1434a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1435a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1436186905e3SBarry Smith       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1437a703fe33SLois Curfman McInnes     }
14389b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
143929bbc08cSBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
144029bbc08cSBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
1441a8c6a408SBarry Smith     if (!snes->computeumfunction) {
144229bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetMinimizationFunction() first");
1443a8c6a408SBarry Smith     }
144429bbc08cSBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetHessian()");
1445d4bb536fSBarry Smith   } else {
144629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
1447d4bb536fSBarry Smith   }
1448a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
144982bf6240SBarry Smith   snes->setupcalled = 1;
14503a40ed3dSBarry Smith   PetscFunctionReturn(0);
14519b94acceSBarry Smith }
14529b94acceSBarry Smith 
14534a2ae208SSatish Balay #undef __FUNCT__
14544a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
14559b94acceSBarry Smith /*@C
14569b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
14579b94acceSBarry Smith    with SNESCreate().
14589b94acceSBarry Smith 
1459c7afd0dbSLois Curfman McInnes    Collective on SNES
1460c7afd0dbSLois Curfman McInnes 
14619b94acceSBarry Smith    Input Parameter:
14629b94acceSBarry Smith .  snes - the SNES context
14639b94acceSBarry Smith 
146436851e7fSLois Curfman McInnes    Level: beginner
146536851e7fSLois Curfman McInnes 
14669b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
14679b94acceSBarry Smith 
146863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
14699b94acceSBarry Smith @*/
14709b94acceSBarry Smith int SNESDestroy(SNES snes)
14719b94acceSBarry Smith {
1472b8a78c4aSBarry Smith   int i,ierr;
14733a40ed3dSBarry Smith 
14743a40ed3dSBarry Smith   PetscFunctionBegin;
147577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
14763a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1477d4bb536fSBarry Smith 
1478be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
14790f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1480be0abb6dSBarry Smith 
1481e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1482606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
14833a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
14843a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
14859b94acceSBarry Smith   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1486522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1487b8a78c4aSBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1488b8a78c4aSBarry Smith     if (snes->monitordestroy[i]) {
1489b8a78c4aSBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1490b8a78c4aSBarry Smith     }
1491b8a78c4aSBarry Smith   }
1492b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)snes);
14930452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
14943a40ed3dSBarry Smith   PetscFunctionReturn(0);
14959b94acceSBarry Smith }
14969b94acceSBarry Smith 
14979b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
14989b94acceSBarry Smith 
14994a2ae208SSatish Balay #undef __FUNCT__
15004a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
15019b94acceSBarry Smith /*@
1502d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
15039b94acceSBarry Smith 
1504c7afd0dbSLois Curfman McInnes    Collective on SNES
1505c7afd0dbSLois Curfman McInnes 
15069b94acceSBarry Smith    Input Parameters:
1507c7afd0dbSLois Curfman McInnes +  snes - the SNES context
150833174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
150933174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
151033174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
151133174efeSLois Curfman McInnes            of the change in the solution between steps
151233174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1513c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1514fee21e36SBarry Smith 
151533174efeSLois Curfman McInnes    Options Database Keys:
1516c7afd0dbSLois Curfman McInnes +    -snes_atol <atol> - Sets atol
1517c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1518c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1519c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1520c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
15219b94acceSBarry Smith 
1522d7a720efSLois Curfman McInnes    Notes:
15239b94acceSBarry Smith    The default maximum number of iterations is 50.
15249b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
15259b94acceSBarry Smith 
152636851e7fSLois Curfman McInnes    Level: intermediate
152736851e7fSLois Curfman McInnes 
152833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
15299b94acceSBarry Smith 
1530d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
15319b94acceSBarry Smith @*/
1532329f5518SBarry Smith int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
15339b94acceSBarry Smith {
15343a40ed3dSBarry Smith   PetscFunctionBegin;
153577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1536d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1537d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1538d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1539d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1540d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
15413a40ed3dSBarry Smith   PetscFunctionReturn(0);
15429b94acceSBarry Smith }
15439b94acceSBarry Smith 
15444a2ae208SSatish Balay #undef __FUNCT__
15454a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
15469b94acceSBarry Smith /*@
154733174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
154833174efeSLois Curfman McInnes 
1549c7afd0dbSLois Curfman McInnes    Not Collective
1550c7afd0dbSLois Curfman McInnes 
155133174efeSLois Curfman McInnes    Input Parameters:
1552c7afd0dbSLois Curfman McInnes +  snes - the SNES context
155333174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
155433174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
155533174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
155633174efeSLois Curfman McInnes            of the change in the solution between steps
155733174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1558c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1559fee21e36SBarry Smith 
156033174efeSLois Curfman McInnes    Notes:
156133174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
156233174efeSLois Curfman McInnes 
156336851e7fSLois Curfman McInnes    Level: intermediate
156436851e7fSLois Curfman McInnes 
156533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
156633174efeSLois Curfman McInnes 
156733174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
156833174efeSLois Curfman McInnes @*/
1569329f5518SBarry Smith int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
157033174efeSLois Curfman McInnes {
15713a40ed3dSBarry Smith   PetscFunctionBegin;
157233174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
157333174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
157433174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
157533174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
157633174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
157733174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
15783a40ed3dSBarry Smith   PetscFunctionReturn(0);
157933174efeSLois Curfman McInnes }
158033174efeSLois Curfman McInnes 
15814a2ae208SSatish Balay #undef __FUNCT__
15824a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
158333174efeSLois Curfman McInnes /*@
15849b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
15859b94acceSBarry Smith 
1586fee21e36SBarry Smith    Collective on SNES
1587fee21e36SBarry Smith 
1588c7afd0dbSLois Curfman McInnes    Input Parameters:
1589c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1590c7afd0dbSLois Curfman McInnes -  tol - tolerance
1591c7afd0dbSLois Curfman McInnes 
15929b94acceSBarry Smith    Options Database Key:
1593c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
15949b94acceSBarry Smith 
159536851e7fSLois Curfman McInnes    Level: intermediate
159636851e7fSLois Curfman McInnes 
15979b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
15989b94acceSBarry Smith 
1599d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
16009b94acceSBarry Smith @*/
1601329f5518SBarry Smith int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
16029b94acceSBarry Smith {
16033a40ed3dSBarry Smith   PetscFunctionBegin;
160477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16059b94acceSBarry Smith   snes->deltatol = tol;
16063a40ed3dSBarry Smith   PetscFunctionReturn(0);
16079b94acceSBarry Smith }
16089b94acceSBarry Smith 
16094a2ae208SSatish Balay #undef __FUNCT__
16104a2ae208SSatish Balay #define __FUNCT__ "SNESSetMinimizationFunctionTolerance"
16119b94acceSBarry Smith /*@
161277c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
16139b94acceSBarry Smith    for unconstrained minimization solvers.
16149b94acceSBarry Smith 
1615fee21e36SBarry Smith    Collective on SNES
1616fee21e36SBarry Smith 
1617c7afd0dbSLois Curfman McInnes    Input Parameters:
1618c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1619c7afd0dbSLois Curfman McInnes -  ftol - minimum function tolerance
1620c7afd0dbSLois Curfman McInnes 
16219b94acceSBarry Smith    Options Database Key:
1622c7afd0dbSLois Curfman McInnes .  -snes_fmin <ftol> - Sets ftol
16239b94acceSBarry Smith 
16249b94acceSBarry Smith    Note:
162577c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
16269b94acceSBarry Smith    methods only.
16279b94acceSBarry Smith 
162836851e7fSLois Curfman McInnes    Level: intermediate
162936851e7fSLois Curfman McInnes 
16309b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
16319b94acceSBarry Smith 
1632d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
16339b94acceSBarry Smith @*/
1634329f5518SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,PetscReal ftol)
16359b94acceSBarry Smith {
16363a40ed3dSBarry Smith   PetscFunctionBegin;
163777c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16389b94acceSBarry Smith   snes->fmin = ftol;
16393a40ed3dSBarry Smith   PetscFunctionReturn(0);
16409b94acceSBarry Smith }
1641df9fa365SBarry Smith /*
1642df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1643df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1644df9fa365SBarry Smith    macros instead of functions
1645df9fa365SBarry Smith */
16464a2ae208SSatish Balay #undef __FUNCT__
16474a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
1648329f5518SBarry Smith int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1649ce1608b8SBarry Smith {
1650ce1608b8SBarry Smith   int ierr;
1651ce1608b8SBarry Smith 
1652ce1608b8SBarry Smith   PetscFunctionBegin;
1653184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1654ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1655ce1608b8SBarry Smith   PetscFunctionReturn(0);
1656ce1608b8SBarry Smith }
1657ce1608b8SBarry Smith 
16584a2ae208SSatish Balay #undef __FUNCT__
16594a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
1660b0a32e0cSBarry Smith int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw)
1661df9fa365SBarry Smith {
1662df9fa365SBarry Smith   int ierr;
1663df9fa365SBarry Smith 
1664df9fa365SBarry Smith   PetscFunctionBegin;
1665df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1666df9fa365SBarry Smith   PetscFunctionReturn(0);
1667df9fa365SBarry Smith }
1668df9fa365SBarry Smith 
16694a2ae208SSatish Balay #undef __FUNCT__
16704a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
1671b0a32e0cSBarry Smith int SNESLGMonitorDestroy(PetscDrawLG draw)
1672df9fa365SBarry Smith {
1673df9fa365SBarry Smith   int ierr;
1674df9fa365SBarry Smith 
1675df9fa365SBarry Smith   PetscFunctionBegin;
1676df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1677df9fa365SBarry Smith   PetscFunctionReturn(0);
1678df9fa365SBarry Smith }
1679df9fa365SBarry Smith 
16809b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
16819b94acceSBarry Smith 
16824a2ae208SSatish Balay #undef __FUNCT__
16834a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
16849b94acceSBarry Smith /*@C
16855cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
16869b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
16879b94acceSBarry Smith    progress.
16889b94acceSBarry Smith 
1689fee21e36SBarry Smith    Collective on SNES
1690fee21e36SBarry Smith 
1691c7afd0dbSLois Curfman McInnes    Input Parameters:
1692c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1693c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1694b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1695b3006f0bSLois Curfman McInnes           monitor routine (use PETSC_NULL if no context is desitre)
1696b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1697b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
16989b94acceSBarry Smith 
1699c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1700329f5518SBarry Smith $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1701c7afd0dbSLois Curfman McInnes 
1702c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1703c7afd0dbSLois Curfman McInnes .    its - iteration number
1704c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
1705c7afd0dbSLois Curfman McInnes             for SNES_NONLINEAR_EQUATIONS methods
170640a0c1c6SLois Curfman McInnes .    norm - 2-norm gradient value (may be estimated)
1707c7afd0dbSLois Curfman McInnes             for SNES_UNCONSTRAINED_MINIMIZATION methods
170840a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
17099b94acceSBarry Smith 
17109665c990SLois Curfman McInnes    Options Database Keys:
1711c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1712c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1713c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1714c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1715c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1716c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1717c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1718c7afd0dbSLois Curfman McInnes                             the options database.
17199665c990SLois Curfman McInnes 
1720639f9d9dSBarry Smith    Notes:
17216bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
17226bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
17236bc08f3fSLois Curfman McInnes    order in which they were set.
1724639f9d9dSBarry Smith 
172536851e7fSLois Curfman McInnes    Level: intermediate
172636851e7fSLois Curfman McInnes 
17279b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
17289b94acceSBarry Smith 
17295cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
17309b94acceSBarry Smith @*/
1731329f5518SBarry Smith int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
17329b94acceSBarry Smith {
17333a40ed3dSBarry Smith   PetscFunctionBegin;
1734184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1735639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
173629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1737639f9d9dSBarry Smith   }
1738639f9d9dSBarry Smith 
1739639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1740b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1741639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
17423a40ed3dSBarry Smith   PetscFunctionReturn(0);
17439b94acceSBarry Smith }
17449b94acceSBarry Smith 
17454a2ae208SSatish Balay #undef __FUNCT__
17464a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
17475cd90555SBarry Smith /*@C
17485cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
17495cd90555SBarry Smith 
1750c7afd0dbSLois Curfman McInnes    Collective on SNES
1751c7afd0dbSLois Curfman McInnes 
17525cd90555SBarry Smith    Input Parameters:
17535cd90555SBarry Smith .  snes - the SNES context
17545cd90555SBarry Smith 
17555cd90555SBarry Smith    Options Database:
1756c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1757c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1758c7afd0dbSLois Curfman McInnes     set via the options database
17595cd90555SBarry Smith 
17605cd90555SBarry Smith    Notes:
17615cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
17625cd90555SBarry Smith 
176336851e7fSLois Curfman McInnes    Level: intermediate
176436851e7fSLois Curfman McInnes 
17655cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
17665cd90555SBarry Smith 
17675cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
17685cd90555SBarry Smith @*/
17695cd90555SBarry Smith int SNESClearMonitor(SNES snes)
17705cd90555SBarry Smith {
17715cd90555SBarry Smith   PetscFunctionBegin;
1772184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17735cd90555SBarry Smith   snes->numbermonitors = 0;
17745cd90555SBarry Smith   PetscFunctionReturn(0);
17755cd90555SBarry Smith }
17765cd90555SBarry Smith 
17774a2ae208SSatish Balay #undef __FUNCT__
17784a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
17799b94acceSBarry Smith /*@C
17809b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
17819b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
17829b94acceSBarry Smith 
1783fee21e36SBarry Smith    Collective on SNES
1784fee21e36SBarry Smith 
1785c7afd0dbSLois Curfman McInnes    Input Parameters:
1786c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1787c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1788c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1789c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
17909b94acceSBarry Smith 
1791c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1792329f5518SBarry Smith $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1793c7afd0dbSLois Curfman McInnes 
1794c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1795c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1796184914b5SBarry Smith .    reason - reason for convergence/divergence
1797c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
1798c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1799c7afd0dbSLois Curfman McInnes .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1800c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1801c7afd0dbSLois Curfman McInnes -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
18029b94acceSBarry Smith 
180336851e7fSLois Curfman McInnes    Level: advanced
180436851e7fSLois Curfman McInnes 
18059b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
18069b94acceSBarry Smith 
180740191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
180840191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
18099b94acceSBarry Smith @*/
1810329f5518SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
18119b94acceSBarry Smith {
18123a40ed3dSBarry Smith   PetscFunctionBegin;
1813184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
18149b94acceSBarry Smith   (snes)->converged = func;
18159b94acceSBarry Smith   (snes)->cnvP      = cctx;
18163a40ed3dSBarry Smith   PetscFunctionReturn(0);
18179b94acceSBarry Smith }
18189b94acceSBarry Smith 
18194a2ae208SSatish Balay #undef __FUNCT__
18204a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
1821184914b5SBarry Smith /*@C
1822184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1823184914b5SBarry Smith 
1824184914b5SBarry Smith    Not Collective
1825184914b5SBarry Smith 
1826184914b5SBarry Smith    Input Parameter:
1827184914b5SBarry Smith .  snes - the SNES context
1828184914b5SBarry Smith 
1829184914b5SBarry Smith    Output Parameter:
1830e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1831184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1832184914b5SBarry Smith 
1833184914b5SBarry Smith    Level: intermediate
1834184914b5SBarry Smith 
1835184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1836184914b5SBarry Smith 
1837184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1838184914b5SBarry Smith 
1839184914b5SBarry Smith .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1840435da068SBarry Smith           SNESConverged_UM_LS(), SNESConverged_UM_TR(), SNESConvergedReason
1841184914b5SBarry Smith @*/
1842184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1843184914b5SBarry Smith {
1844184914b5SBarry Smith   PetscFunctionBegin;
1845184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1846184914b5SBarry Smith   *reason = snes->reason;
1847184914b5SBarry Smith   PetscFunctionReturn(0);
1848184914b5SBarry Smith }
1849184914b5SBarry Smith 
18504a2ae208SSatish Balay #undef __FUNCT__
18514a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1852c9005455SLois Curfman McInnes /*@
1853c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1854c9005455SLois Curfman McInnes 
1855fee21e36SBarry Smith    Collective on SNES
1856fee21e36SBarry Smith 
1857c7afd0dbSLois Curfman McInnes    Input Parameters:
1858c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1859c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1860cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1861758f92a0SBarry Smith .  na  - size of a and its
186264731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1863758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1864c7afd0dbSLois Curfman McInnes 
1865c9005455SLois Curfman McInnes    Notes:
1866c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1867c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1868c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1869c9005455SLois Curfman McInnes    at each step.
1870c9005455SLois Curfman McInnes 
1871c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1872c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1873c9005455SLois Curfman McInnes    during the section of code that is being timed.
1874c9005455SLois Curfman McInnes 
187536851e7fSLois Curfman McInnes    Level: intermediate
187636851e7fSLois Curfman McInnes 
1877c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1878758f92a0SBarry Smith 
187908405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1880758f92a0SBarry Smith 
1881c9005455SLois Curfman McInnes @*/
1882329f5518SBarry Smith int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset)
1883c9005455SLois Curfman McInnes {
18843a40ed3dSBarry Smith   PetscFunctionBegin;
1885c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1886c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1887c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1888758f92a0SBarry Smith   snes->conv_hist_its   = its;
1889758f92a0SBarry Smith   snes->conv_hist_max   = na;
1890758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1891758f92a0SBarry Smith   PetscFunctionReturn(0);
1892758f92a0SBarry Smith }
1893758f92a0SBarry Smith 
18944a2ae208SSatish Balay #undef __FUNCT__
18954a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
18960c4c9dddSBarry Smith /*@C
1897758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1898758f92a0SBarry Smith 
1899758f92a0SBarry Smith    Collective on SNES
1900758f92a0SBarry Smith 
1901758f92a0SBarry Smith    Input Parameter:
1902758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1903758f92a0SBarry Smith 
1904758f92a0SBarry Smith    Output Parameters:
1905758f92a0SBarry Smith .  a   - array to hold history
1906758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1907758f92a0SBarry Smith          negative if not converged) for each solve.
1908758f92a0SBarry Smith -  na  - size of a and its
1909758f92a0SBarry Smith 
1910758f92a0SBarry Smith    Notes:
1911758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1912758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1913758f92a0SBarry Smith 
1914758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1915758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1916758f92a0SBarry Smith    during the section of code that is being timed.
1917758f92a0SBarry Smith 
1918758f92a0SBarry Smith    Level: intermediate
1919758f92a0SBarry Smith 
1920758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1921758f92a0SBarry Smith 
1922758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1923758f92a0SBarry Smith 
1924758f92a0SBarry Smith @*/
1925329f5518SBarry Smith int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na)
1926758f92a0SBarry Smith {
1927758f92a0SBarry Smith   PetscFunctionBegin;
1928758f92a0SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1929758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1930758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1931758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
19323a40ed3dSBarry Smith   PetscFunctionReturn(0);
1933c9005455SLois Curfman McInnes }
1934c9005455SLois Curfman McInnes 
1935e74ef692SMatthew Knepley #undef __FUNCT__
1936e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC"
193776b2cf59SMatthew Knepley /*@
193876b2cf59SMatthew Knepley   SNESSetRhsBC - Sets the function which applies boundary conditions
193976b2cf59SMatthew Knepley   to the Rhs of each system.
194076b2cf59SMatthew Knepley 
194176b2cf59SMatthew Knepley   Collective on SNES
194276b2cf59SMatthew Knepley 
194376b2cf59SMatthew Knepley   Input Parameters:
194476b2cf59SMatthew Knepley . snes - The nonlinear solver context
194576b2cf59SMatthew Knepley . func - The function
194676b2cf59SMatthew Knepley 
194776b2cf59SMatthew Knepley   Calling sequence of func:
194876b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx);
194976b2cf59SMatthew Knepley 
195076b2cf59SMatthew Knepley . rhs - The current rhs vector
195176b2cf59SMatthew Knepley . ctx - The user-context
195276b2cf59SMatthew Knepley 
195376b2cf59SMatthew Knepley   Level: intermediate
195476b2cf59SMatthew Knepley 
195576b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
195676b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
195776b2cf59SMatthew Knepley @*/
195876b2cf59SMatthew Knepley int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *))
195976b2cf59SMatthew Knepley {
196076b2cf59SMatthew Knepley   PetscFunctionBegin;
196176b2cf59SMatthew Knepley   PetscValidHeaderSpecific(snes, SNES_COOKIE);
196276b2cf59SMatthew Knepley   snes->applyrhsbc = func;
196376b2cf59SMatthew Knepley   PetscFunctionReturn(0);
196476b2cf59SMatthew Knepley }
196576b2cf59SMatthew Knepley 
1966e74ef692SMatthew Knepley #undef __FUNCT__
1967e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC"
196876b2cf59SMatthew Knepley /*@
196976b2cf59SMatthew Knepley   SNESDefaultRhsBC - The default boundary condition function which does nothing.
197076b2cf59SMatthew Knepley 
197176b2cf59SMatthew Knepley   Not collective
197276b2cf59SMatthew Knepley 
197376b2cf59SMatthew Knepley   Input Parameters:
197476b2cf59SMatthew Knepley . snes - The nonlinear solver context
197576b2cf59SMatthew Knepley . rhs  - The Rhs
197676b2cf59SMatthew Knepley . ctx  - The user-context
197776b2cf59SMatthew Knepley 
197876b2cf59SMatthew Knepley   Level: beginner
197976b2cf59SMatthew Knepley 
198076b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
198176b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
198276b2cf59SMatthew Knepley @*/
198376b2cf59SMatthew Knepley int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
198476b2cf59SMatthew Knepley {
198576b2cf59SMatthew Knepley   PetscFunctionBegin;
198676b2cf59SMatthew Knepley   PetscFunctionReturn(0);
198776b2cf59SMatthew Knepley }
198876b2cf59SMatthew Knepley 
1989e74ef692SMatthew Knepley #undef __FUNCT__
1990e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC"
199176b2cf59SMatthew Knepley /*@
199276b2cf59SMatthew Knepley   SNESSetSolutionBC - Sets the function which applies boundary conditions
199376b2cf59SMatthew Knepley   to the solution of each system.
199476b2cf59SMatthew Knepley 
199576b2cf59SMatthew Knepley   Collective on SNES
199676b2cf59SMatthew Knepley 
199776b2cf59SMatthew Knepley   Input Parameters:
199876b2cf59SMatthew Knepley . snes - The nonlinear solver context
199976b2cf59SMatthew Knepley . func - The function
200076b2cf59SMatthew Knepley 
200176b2cf59SMatthew Knepley   Calling sequence of func:
200276b2cf59SMatthew Knepley . func (TS ts, Vec rsol, void *ctx);
200376b2cf59SMatthew Knepley 
200476b2cf59SMatthew Knepley . sol - The current solution vector
200576b2cf59SMatthew Knepley . ctx - The user-context
200676b2cf59SMatthew Knepley 
200776b2cf59SMatthew Knepley   Level: intermediate
200876b2cf59SMatthew Knepley 
200976b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
201076b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
201176b2cf59SMatthew Knepley @*/
201276b2cf59SMatthew Knepley int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *))
201376b2cf59SMatthew Knepley {
201476b2cf59SMatthew Knepley   PetscFunctionBegin;
201576b2cf59SMatthew Knepley   PetscValidHeaderSpecific(snes, SNES_COOKIE);
201676b2cf59SMatthew Knepley   snes->applysolbc = func;
201776b2cf59SMatthew Knepley   PetscFunctionReturn(0);
201876b2cf59SMatthew Knepley }
201976b2cf59SMatthew Knepley 
2020e74ef692SMatthew Knepley #undef __FUNCT__
2021e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC"
202276b2cf59SMatthew Knepley /*@
202376b2cf59SMatthew Knepley   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
202476b2cf59SMatthew Knepley 
202576b2cf59SMatthew Knepley   Not collective
202676b2cf59SMatthew Knepley 
202776b2cf59SMatthew Knepley   Input Parameters:
202876b2cf59SMatthew Knepley . snes - The nonlinear solver context
202976b2cf59SMatthew Knepley . sol  - The solution
203076b2cf59SMatthew Knepley . ctx  - The user-context
203176b2cf59SMatthew Knepley 
203276b2cf59SMatthew Knepley   Level: beginner
203376b2cf59SMatthew Knepley 
203476b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
203576b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
203676b2cf59SMatthew Knepley @*/
203776b2cf59SMatthew Knepley int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
203876b2cf59SMatthew Knepley {
203976b2cf59SMatthew Knepley   PetscFunctionBegin;
204076b2cf59SMatthew Knepley   PetscFunctionReturn(0);
204176b2cf59SMatthew Knepley }
204276b2cf59SMatthew Knepley 
2043e74ef692SMatthew Knepley #undef __FUNCT__
2044e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
204576b2cf59SMatthew Knepley /*@
204676b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
204776b2cf59SMatthew Knepley   at the beginning of every step of the iteration.
204876b2cf59SMatthew Knepley 
204976b2cf59SMatthew Knepley   Collective on SNES
205076b2cf59SMatthew Knepley 
205176b2cf59SMatthew Knepley   Input Parameters:
205276b2cf59SMatthew Knepley . snes - The nonlinear solver context
205376b2cf59SMatthew Knepley . func - The function
205476b2cf59SMatthew Knepley 
205576b2cf59SMatthew Knepley   Calling sequence of func:
205676b2cf59SMatthew Knepley . func (TS ts, int step);
205776b2cf59SMatthew Knepley 
205876b2cf59SMatthew Knepley . step - The current step of the iteration
205976b2cf59SMatthew Knepley 
206076b2cf59SMatthew Knepley   Level: intermediate
206176b2cf59SMatthew Knepley 
206276b2cf59SMatthew Knepley .keywords: SNES, update
206376b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
206476b2cf59SMatthew Knepley @*/
206576b2cf59SMatthew Knepley int SNESSetUpdate(SNES snes, int (*func)(SNES, int))
206676b2cf59SMatthew Knepley {
206776b2cf59SMatthew Knepley   PetscFunctionBegin;
206876b2cf59SMatthew Knepley   PetscValidHeaderSpecific(snes, SNES_COOKIE);
206976b2cf59SMatthew Knepley   snes->update = func;
207076b2cf59SMatthew Knepley   PetscFunctionReturn(0);
207176b2cf59SMatthew Knepley }
207276b2cf59SMatthew Knepley 
2073e74ef692SMatthew Knepley #undef __FUNCT__
2074e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
207576b2cf59SMatthew Knepley /*@
207676b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
207776b2cf59SMatthew Knepley 
207876b2cf59SMatthew Knepley   Not collective
207976b2cf59SMatthew Knepley 
208076b2cf59SMatthew Knepley   Input Parameters:
208176b2cf59SMatthew Knepley . snes - The nonlinear solver context
208276b2cf59SMatthew Knepley . step - The current step of the iteration
208376b2cf59SMatthew Knepley 
208476b2cf59SMatthew Knepley .keywords: SNES, update
208576b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
208676b2cf59SMatthew Knepley @*/
208776b2cf59SMatthew Knepley int SNESDefaultUpdate(SNES snes, int step)
208876b2cf59SMatthew Knepley {
208976b2cf59SMatthew Knepley   PetscFunctionBegin;
209076b2cf59SMatthew Knepley   PetscFunctionReturn(0);
209176b2cf59SMatthew Knepley }
209276b2cf59SMatthew Knepley 
20934a2ae208SSatish Balay #undef __FUNCT__
20944a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
20959b94acceSBarry Smith /*
20969b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
20979b94acceSBarry Smith    positive parameter delta.
20989b94acceSBarry Smith 
20999b94acceSBarry Smith     Input Parameters:
2100c7afd0dbSLois Curfman McInnes +   snes - the SNES context
21019b94acceSBarry Smith .   y - approximate solution of linear system
21029b94acceSBarry Smith .   fnorm - 2-norm of current function
2103c7afd0dbSLois Curfman McInnes -   delta - trust region size
21049b94acceSBarry Smith 
21059b94acceSBarry Smith     Output Parameters:
2106c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
21079b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
21089b94acceSBarry Smith     region, and exceeds zero otherwise.
2109c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
21109b94acceSBarry Smith 
21119b94acceSBarry Smith     Note:
21126831982aSBarry Smith     For non-trust region methods such as SNESEQLS, the parameter delta
21139b94acceSBarry Smith     is set to be the maximum allowable step size.
21149b94acceSBarry Smith 
21159b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
21169b94acceSBarry Smith */
2117329f5518SBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,
2118329f5518SBarry Smith                           PetscReal *gpnorm,PetscReal *ynorm)
21199b94acceSBarry Smith {
2120329f5518SBarry Smith   PetscReal norm;
2121ea709b57SSatish Balay   PetscScalar cnorm;
21223a40ed3dSBarry Smith   int    ierr;
21233a40ed3dSBarry Smith 
21243a40ed3dSBarry Smith   PetscFunctionBegin;
2125184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2126184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
2127184914b5SBarry Smith   PetscCheckSameComm(snes,y);
2128184914b5SBarry Smith 
21293a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
21309b94acceSBarry Smith   if (norm > *delta) {
21319b94acceSBarry Smith      norm = *delta/norm;
21329b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
21339b94acceSBarry Smith      cnorm = norm;
2134329f5518SBarry Smith      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
21359b94acceSBarry Smith      *ynorm = *delta;
21369b94acceSBarry Smith   } else {
21379b94acceSBarry Smith      *gpnorm = 0.0;
21389b94acceSBarry Smith      *ynorm = norm;
21399b94acceSBarry Smith   }
21403a40ed3dSBarry Smith   PetscFunctionReturn(0);
21419b94acceSBarry Smith }
21429b94acceSBarry Smith 
21434a2ae208SSatish Balay #undef __FUNCT__
21444a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
21459b94acceSBarry Smith /*@
21469b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
214763a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
21489b94acceSBarry Smith 
2149c7afd0dbSLois Curfman McInnes    Collective on SNES
2150c7afd0dbSLois Curfman McInnes 
2151b2002411SLois Curfman McInnes    Input Parameters:
2152c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2153c7afd0dbSLois Curfman McInnes -  x - the solution vector
21549b94acceSBarry Smith 
21559b94acceSBarry Smith    Output Parameter:
2156b2002411SLois Curfman McInnes .  its - number of iterations until termination
21579b94acceSBarry Smith 
2158b2002411SLois Curfman McInnes    Notes:
21598ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
21608ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
21618ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
21628ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
21638ddd3da0SLois Curfman McInnes 
216436851e7fSLois Curfman McInnes    Level: beginner
216536851e7fSLois Curfman McInnes 
21669b94acceSBarry Smith .keywords: SNES, nonlinear, solve
21679b94acceSBarry Smith 
216863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
21699b94acceSBarry Smith @*/
21708ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
21719b94acceSBarry Smith {
2172f1af5d2fSBarry Smith   int        ierr;
2173f1af5d2fSBarry Smith   PetscTruth flg;
2174052efed2SBarry Smith 
21753a40ed3dSBarry Smith   PetscFunctionBegin;
217677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2177184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
2178184914b5SBarry Smith   PetscCheckSameComm(snes,x);
217974679c65SBarry Smith   PetscValidIntPointer(its);
218029bbc08cSBarry Smith   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
218174637425SBarry Smith 
218282bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
2183c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
2184758f92a0SBarry Smith   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
2185*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2186c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
21879b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
2188*d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
2189b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
219093b98531SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); }
21913a40ed3dSBarry Smith   PetscFunctionReturn(0);
21929b94acceSBarry Smith }
21939b94acceSBarry Smith 
21949b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
21959b94acceSBarry Smith 
21964a2ae208SSatish Balay #undef __FUNCT__
21974a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
219882bf6240SBarry Smith /*@C
21994b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
22009b94acceSBarry Smith 
2201fee21e36SBarry Smith    Collective on SNES
2202fee21e36SBarry Smith 
2203c7afd0dbSLois Curfman McInnes    Input Parameters:
2204c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2205454a90a3SBarry Smith -  type - a known method
2206c7afd0dbSLois Curfman McInnes 
2207c7afd0dbSLois Curfman McInnes    Options Database Key:
2208454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
2209c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
2210ae12b187SLois Curfman McInnes 
22119b94acceSBarry Smith    Notes:
2212e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
22136831982aSBarry Smith +    SNESEQLS - Newton's method with line search
2214c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
22156831982aSBarry Smith .    SNESEQTR - Newton's method with trust region
2216c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
22176831982aSBarry Smith .    SNESUMTR - Newton's method with trust region
2218c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
22196831982aSBarry Smith -    SNESUMLS - Newton's method with line search
2220c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
22219b94acceSBarry Smith 
2222ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
2223ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
2224ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
2225ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
2226ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
2227ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
2228ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
2229ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
2230ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
2231b0a32e0cSBarry Smith   appropriate method.
223236851e7fSLois Curfman McInnes 
223336851e7fSLois Curfman McInnes   Level: intermediate
2234a703fe33SLois Curfman McInnes 
2235454a90a3SBarry Smith .keywords: SNES, set, type
2236435da068SBarry Smith 
2237435da068SBarry Smith .seealso: SNESType, SNESCreate()
2238435da068SBarry Smith 
22399b94acceSBarry Smith @*/
2240454a90a3SBarry Smith int SNESSetType(SNES snes,SNESType type)
22419b94acceSBarry Smith {
22420f5bd95cSBarry Smith   int        ierr,(*r)(SNES);
22436831982aSBarry Smith   PetscTruth match;
22443a40ed3dSBarry Smith 
22453a40ed3dSBarry Smith   PetscFunctionBegin;
224677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
22470f5bd95cSBarry Smith   PetscValidCharPointer(type);
224882bf6240SBarry Smith 
22496831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
22500f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
225182bf6240SBarry Smith 
225282bf6240SBarry Smith   if (snes->setupcalled) {
2253e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
225482bf6240SBarry Smith     snes->data = 0;
225582bf6240SBarry Smith   }
22567d1a2b2bSBarry Smith 
22579b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
225882bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
225982bf6240SBarry Smith 
2260b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
226182bf6240SBarry Smith 
226229bbc08cSBarry Smith   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
22631548b14aSBarry Smith 
2264606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
226582bf6240SBarry Smith   snes->data = 0;
22663a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
226782bf6240SBarry Smith 
2268454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
226982bf6240SBarry Smith   snes->set_method_called = 1;
227082bf6240SBarry Smith 
22713a40ed3dSBarry Smith   PetscFunctionReturn(0);
22729b94acceSBarry Smith }
22739b94acceSBarry Smith 
2274a847f771SSatish Balay 
22759b94acceSBarry Smith /* --------------------------------------------------------------------- */
22764a2ae208SSatish Balay #undef __FUNCT__
22774a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
22789b94acceSBarry Smith /*@C
22799b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2280f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
22819b94acceSBarry Smith 
2282fee21e36SBarry Smith    Not Collective
2283fee21e36SBarry Smith 
228436851e7fSLois Curfman McInnes    Level: advanced
228536851e7fSLois Curfman McInnes 
22869b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
22879b94acceSBarry Smith 
22889b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
22899b94acceSBarry Smith @*/
2290cf256101SBarry Smith int SNESRegisterDestroy(void)
22919b94acceSBarry Smith {
229282bf6240SBarry Smith   int ierr;
229382bf6240SBarry Smith 
22943a40ed3dSBarry Smith   PetscFunctionBegin;
229582bf6240SBarry Smith   if (SNESList) {
2296b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
229782bf6240SBarry Smith     SNESList = 0;
22989b94acceSBarry Smith   }
22994c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
23003a40ed3dSBarry Smith   PetscFunctionReturn(0);
23019b94acceSBarry Smith }
23029b94acceSBarry Smith 
23034a2ae208SSatish Balay #undef __FUNCT__
23044a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
23059b94acceSBarry Smith /*@C
23069a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
23079b94acceSBarry Smith 
2308c7afd0dbSLois Curfman McInnes    Not Collective
2309c7afd0dbSLois Curfman McInnes 
23109b94acceSBarry Smith    Input Parameter:
23114b0e389bSBarry Smith .  snes - nonlinear solver context
23129b94acceSBarry Smith 
23139b94acceSBarry Smith    Output Parameter:
23143a7fca6bSBarry Smith .  type - SNES method (a character string)
23159b94acceSBarry Smith 
231636851e7fSLois Curfman McInnes    Level: intermediate
231736851e7fSLois Curfman McInnes 
2318454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
23199b94acceSBarry Smith @*/
2320454a90a3SBarry Smith int SNESGetType(SNES snes,SNESType *type)
23219b94acceSBarry Smith {
23223a40ed3dSBarry Smith   PetscFunctionBegin;
2323184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2324454a90a3SBarry Smith   *type = snes->type_name;
23253a40ed3dSBarry Smith   PetscFunctionReturn(0);
23269b94acceSBarry Smith }
23279b94acceSBarry Smith 
23284a2ae208SSatish Balay #undef __FUNCT__
23294a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
23309b94acceSBarry Smith /*@C
23319b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
23329b94acceSBarry Smith    stored.
23339b94acceSBarry Smith 
2334c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2335c7afd0dbSLois Curfman McInnes 
23369b94acceSBarry Smith    Input Parameter:
23379b94acceSBarry Smith .  snes - the SNES context
23389b94acceSBarry Smith 
23399b94acceSBarry Smith    Output Parameter:
23409b94acceSBarry Smith .  x - the solution
23419b94acceSBarry Smith 
234236851e7fSLois Curfman McInnes    Level: advanced
234336851e7fSLois Curfman McInnes 
23449b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
23459b94acceSBarry Smith 
2346a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
23479b94acceSBarry Smith @*/
23489b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
23499b94acceSBarry Smith {
23503a40ed3dSBarry Smith   PetscFunctionBegin;
235177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
23529b94acceSBarry Smith   *x = snes->vec_sol_always;
23533a40ed3dSBarry Smith   PetscFunctionReturn(0);
23549b94acceSBarry Smith }
23559b94acceSBarry Smith 
23564a2ae208SSatish Balay #undef __FUNCT__
23574a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
23589b94acceSBarry Smith /*@C
23599b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
23609b94acceSBarry Smith    stored.
23619b94acceSBarry Smith 
2362c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2363c7afd0dbSLois Curfman McInnes 
23649b94acceSBarry Smith    Input Parameter:
23659b94acceSBarry Smith .  snes - the SNES context
23669b94acceSBarry Smith 
23679b94acceSBarry Smith    Output Parameter:
23689b94acceSBarry Smith .  x - the solution update
23699b94acceSBarry Smith 
237036851e7fSLois Curfman McInnes    Level: advanced
237136851e7fSLois Curfman McInnes 
23729b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
23739b94acceSBarry Smith 
23749b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
23759b94acceSBarry Smith @*/
23769b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
23779b94acceSBarry Smith {
23783a40ed3dSBarry Smith   PetscFunctionBegin;
237977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
23809b94acceSBarry Smith   *x = snes->vec_sol_update_always;
23813a40ed3dSBarry Smith   PetscFunctionReturn(0);
23829b94acceSBarry Smith }
23839b94acceSBarry Smith 
23844a2ae208SSatish Balay #undef __FUNCT__
23854a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
23869b94acceSBarry Smith /*@C
23873638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
23889b94acceSBarry Smith 
2389c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2390c7afd0dbSLois Curfman McInnes 
23919b94acceSBarry Smith    Input Parameter:
23929b94acceSBarry Smith .  snes - the SNES context
23939b94acceSBarry Smith 
23949b94acceSBarry Smith    Output Parameter:
23957bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
239600036973SBarry Smith .  ctx - the function context (or PETSC_NULL)
239700036973SBarry Smith -  func - the function (or PETSC_NULL)
23989b94acceSBarry Smith 
23999b94acceSBarry Smith    Notes:
24009b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
24019b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
24029b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
24039b94acceSBarry Smith 
240436851e7fSLois Curfman McInnes    Level: advanced
240536851e7fSLois Curfman McInnes 
2406a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
24079b94acceSBarry Smith 
240831615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
240931615d04SLois Curfman McInnes           SNESGetGradient()
24107bf4e008SBarry Smith 
24119b94acceSBarry Smith @*/
241200036973SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
24139b94acceSBarry Smith {
24143a40ed3dSBarry Smith   PetscFunctionBegin;
241577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2416a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
241729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
2418a8c6a408SBarry Smith   }
24197bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
24207bf4e008SBarry Smith   if (ctx)  *ctx  = snes->funP;
242100036973SBarry Smith   if (func) *func = snes->computefunction;
24223a40ed3dSBarry Smith   PetscFunctionReturn(0);
24239b94acceSBarry Smith }
24249b94acceSBarry Smith 
24254a2ae208SSatish Balay #undef __FUNCT__
24264a2ae208SSatish Balay #define __FUNCT__ "SNESGetGradient"
24279b94acceSBarry Smith /*@C
24283638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
24299b94acceSBarry Smith 
2430c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2431c7afd0dbSLois Curfman McInnes 
24329b94acceSBarry Smith    Input Parameter:
24339b94acceSBarry Smith .  snes - the SNES context
24349b94acceSBarry Smith 
24359b94acceSBarry Smith    Output Parameter:
24367bf4e008SBarry Smith +  r - the gradient (or PETSC_NULL)
24377bf4e008SBarry Smith -  ctx - the gradient context (or PETSC_NULL)
24389b94acceSBarry Smith 
24399b94acceSBarry Smith    Notes:
24409b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
24419b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
24429b94acceSBarry Smith    SNESGetFunction().
24439b94acceSBarry Smith 
244436851e7fSLois Curfman McInnes    Level: advanced
244536851e7fSLois Curfman McInnes 
24469b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
24479b94acceSBarry Smith 
24487bf4e008SBarry Smith .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(),
24497bf4e008SBarry Smith           SNESSetGradient(), SNESSetFunction()
24507bf4e008SBarry Smith 
24519b94acceSBarry Smith @*/
24527bf4e008SBarry Smith int SNESGetGradient(SNES snes,Vec *r,void **ctx)
24539b94acceSBarry Smith {
24543a40ed3dSBarry Smith   PetscFunctionBegin;
245577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
24563a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
245729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
24583a40ed3dSBarry Smith   }
24597bf4e008SBarry Smith   if (r)   *r = snes->vec_func_always;
24607bf4e008SBarry Smith   if (ctx) *ctx = snes->funP;
24613a40ed3dSBarry Smith   PetscFunctionReturn(0);
24629b94acceSBarry Smith }
24639b94acceSBarry Smith 
24644a2ae208SSatish Balay #undef __FUNCT__
24654a2ae208SSatish Balay #define __FUNCT__ "SNESGetMinimizationFunction"
24667bf4e008SBarry Smith /*@C
24679b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
24689b94acceSBarry Smith    unconstrained minimization problems.
24699b94acceSBarry Smith 
2470c7afd0dbSLois Curfman McInnes    Not Collective
2471c7afd0dbSLois Curfman McInnes 
24729b94acceSBarry Smith    Input Parameter:
24739b94acceSBarry Smith .  snes - the SNES context
24749b94acceSBarry Smith 
24759b94acceSBarry Smith    Output Parameter:
24767bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
24777bf4e008SBarry Smith -  ctx - the function context (or PETSC_NULL)
24789b94acceSBarry Smith 
24799b94acceSBarry Smith    Notes:
24809b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
24819b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
24829b94acceSBarry Smith    SNESGetFunction().
24839b94acceSBarry Smith 
248436851e7fSLois Curfman McInnes    Level: advanced
248536851e7fSLois Curfman McInnes 
24869b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
24879b94acceSBarry Smith 
24887bf4e008SBarry Smith .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction()
24897bf4e008SBarry Smith 
24909b94acceSBarry Smith @*/
2491329f5518SBarry Smith int SNESGetMinimizationFunction(SNES snes,PetscReal *r,void **ctx)
24929b94acceSBarry Smith {
24933a40ed3dSBarry Smith   PetscFunctionBegin;
249477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
249574679c65SBarry Smith   PetscValidScalarPointer(r);
24963a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
249729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
24983a40ed3dSBarry Smith   }
24997bf4e008SBarry Smith   if (r)   *r = snes->fc;
25007bf4e008SBarry Smith   if (ctx) *ctx = snes->umfunP;
25013a40ed3dSBarry Smith   PetscFunctionReturn(0);
25029b94acceSBarry Smith }
25039b94acceSBarry Smith 
25044a2ae208SSatish Balay #undef __FUNCT__
25054a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
25063c7409f5SSatish Balay /*@C
25073c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2508d850072dSLois Curfman McInnes    SNES options in the database.
25093c7409f5SSatish Balay 
2510fee21e36SBarry Smith    Collective on SNES
2511fee21e36SBarry Smith 
2512c7afd0dbSLois Curfman McInnes    Input Parameter:
2513c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2514c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2515c7afd0dbSLois Curfman McInnes 
2516d850072dSLois Curfman McInnes    Notes:
2517a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2518c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2519d850072dSLois Curfman McInnes 
252036851e7fSLois Curfman McInnes    Level: advanced
252136851e7fSLois Curfman McInnes 
25223c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2523a86d99e1SLois Curfman McInnes 
2524a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
25253c7409f5SSatish Balay @*/
25263c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
25273c7409f5SSatish Balay {
25283c7409f5SSatish Balay   int ierr;
25293c7409f5SSatish Balay 
25303a40ed3dSBarry Smith   PetscFunctionBegin;
253177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2532639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
25333c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
25343a40ed3dSBarry Smith   PetscFunctionReturn(0);
25353c7409f5SSatish Balay }
25363c7409f5SSatish Balay 
25374a2ae208SSatish Balay #undef __FUNCT__
25384a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
25393c7409f5SSatish Balay /*@C
2540f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2541d850072dSLois Curfman McInnes    SNES options in the database.
25423c7409f5SSatish Balay 
2543fee21e36SBarry Smith    Collective on SNES
2544fee21e36SBarry Smith 
2545c7afd0dbSLois Curfman McInnes    Input Parameters:
2546c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2547c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2548c7afd0dbSLois Curfman McInnes 
2549d850072dSLois Curfman McInnes    Notes:
2550a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2551c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2552d850072dSLois Curfman McInnes 
255336851e7fSLois Curfman McInnes    Level: advanced
255436851e7fSLois Curfman McInnes 
25553c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2556a86d99e1SLois Curfman McInnes 
2557a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
25583c7409f5SSatish Balay @*/
25593c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
25603c7409f5SSatish Balay {
25613c7409f5SSatish Balay   int ierr;
25623c7409f5SSatish Balay 
25633a40ed3dSBarry Smith   PetscFunctionBegin;
256477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2565639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
25663c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
25673a40ed3dSBarry Smith   PetscFunctionReturn(0);
25683c7409f5SSatish Balay }
25693c7409f5SSatish Balay 
25704a2ae208SSatish Balay #undef __FUNCT__
25714a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
25729ab63eb5SSatish Balay /*@C
25733c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
25743c7409f5SSatish Balay    SNES options in the database.
25753c7409f5SSatish Balay 
2576c7afd0dbSLois Curfman McInnes    Not Collective
2577c7afd0dbSLois Curfman McInnes 
25783c7409f5SSatish Balay    Input Parameter:
25793c7409f5SSatish Balay .  snes - the SNES context
25803c7409f5SSatish Balay 
25813c7409f5SSatish Balay    Output Parameter:
25823c7409f5SSatish Balay .  prefix - pointer to the prefix string used
25833c7409f5SSatish Balay 
25849ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
25859ab63eb5SSatish Balay    sufficient length to hold the prefix.
25869ab63eb5SSatish Balay 
258736851e7fSLois Curfman McInnes    Level: advanced
258836851e7fSLois Curfman McInnes 
25893c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2590a86d99e1SLois Curfman McInnes 
2591a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
25923c7409f5SSatish Balay @*/
25933c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
25943c7409f5SSatish Balay {
25953c7409f5SSatish Balay   int ierr;
25963c7409f5SSatish Balay 
25973a40ed3dSBarry Smith   PetscFunctionBegin;
259877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2599639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
26003a40ed3dSBarry Smith   PetscFunctionReturn(0);
26013c7409f5SSatish Balay }
26023c7409f5SSatish Balay 
2603acb85ae6SSatish Balay /*MC
2604f1af5d2fSBarry Smith    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
26059b94acceSBarry Smith 
2606b2002411SLois Curfman McInnes    Synopsis:
26073a7fca6bSBarry Smith    int SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
26089b94acceSBarry Smith 
26098d76a1e5SLois Curfman McInnes    Not collective
26108d76a1e5SLois Curfman McInnes 
2611b2002411SLois Curfman McInnes    Input Parameters:
2612c7afd0dbSLois Curfman McInnes +  name_solver - name of a new user-defined solver
2613b2002411SLois Curfman McInnes .  path - path (either absolute or relative) the library containing this solver
2614b2002411SLois Curfman McInnes .  name_create - name of routine to create method context
2615c7afd0dbSLois Curfman McInnes -  routine_create - routine to create method context
26169b94acceSBarry Smith 
2617b2002411SLois Curfman McInnes    Notes:
2618f1af5d2fSBarry Smith    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
26193a40ed3dSBarry Smith 
2620b2002411SLois Curfman McInnes    If dynamic libraries are used, then the fourth input argument (routine_create)
2621b2002411SLois Curfman McInnes    is ignored.
2622b2002411SLois Curfman McInnes 
2623b9617806SBarry Smith    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT},
26245ba88a07SLois Curfman McInnes    and others of the form ${any_environmental_variable} occuring in pathname will be
26255ba88a07SLois Curfman McInnes    replaced with appropriate values.
26265ba88a07SLois Curfman McInnes 
2627b2002411SLois Curfman McInnes    Sample usage:
262869225d78SLois Curfman McInnes .vb
2629f1af5d2fSBarry Smith    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2630b2002411SLois Curfman McInnes                 "MySolverCreate",MySolverCreate);
263169225d78SLois Curfman McInnes .ve
2632b2002411SLois Curfman McInnes 
2633b2002411SLois Curfman McInnes    Then, your solver can be chosen with the procedural interface via
2634b2002411SLois Curfman McInnes $     SNESSetType(snes,"my_solver")
2635b2002411SLois Curfman McInnes    or at runtime via the option
2636b2002411SLois Curfman McInnes $     -snes_type my_solver
2637b2002411SLois Curfman McInnes 
263836851e7fSLois Curfman McInnes    Level: advanced
263936851e7fSLois Curfman McInnes 
2640b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register
2641b2002411SLois Curfman McInnes 
2642b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2643b2002411SLois Curfman McInnes M*/
2644b2002411SLois Curfman McInnes 
26454a2ae208SSatish Balay #undef __FUNCT__
26464a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
2647f1af5d2fSBarry Smith int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES))
2648b2002411SLois Curfman McInnes {
2649b2002411SLois Curfman McInnes   char fullname[256];
2650b2002411SLois Curfman McInnes   int  ierr;
2651b2002411SLois Curfman McInnes 
2652b2002411SLois Curfman McInnes   PetscFunctionBegin;
2653b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2654c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2655b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2656b2002411SLois Curfman McInnes }
2657